CMS 3D CMS Logo

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

#include <RPCIntegrator.h>

Public Member Functions

L1Phase2MuDTPhDigi createL1Phase2MuDTPhDigi (RPCDetId rpcDetId, int rpc_bx, double rpc_time, double rpc_global_phi, double phiB, int rpc_flag)
 
void finish ()
 
bool hasPosRF_rpc (int wh, int sec) const
 
void initialise (const edm::EventSetup &iEventSetup, double shift_back_fromDT)
 
void makeRPCOnlySegments ()
 
RPCMetaprimitivematchDTwithRPC (cmsdt::metaPrimitive *dt_metaprimitive)
 
void matchWithDTAndUseRPCTime (std::vector< cmsdt::metaPrimitive > &dt_metaprimitives)
 
double phi_DT_MP_conv (double rpc_global_phi, int rpcSector)
 
double phiBending (RPCMetaprimitive *rpc_hit_1, RPCMetaprimitive *rpc_hit_2)
 
int phiInDTTPFormat (double rpc_global_phi, int rpcSector)
 
void prepareMetaPrimitives (edm::Handle< RPCRecHitCollection > rpcRecHits)
 
void removeRPCHitsUsed ()
 
GlobalPoint RPCGlobalPosition (RPCDetId rpcId, const RPCRecHit &rpcIt) const
 
 RPCIntegrator (const edm::ParameterSet &pset, edm::ConsumesCollector &iC)
 
void storeRPCSingleHits ()
 
 ~RPCIntegrator ()
 

Public Attributes

std::vector< RPCMetaprimitiveRPCMetaprimitives_
 
std::vector< L1Phase2MuDTPhDigirpcRecHits_translated_
 

Private Attributes

DTGeometry const * dtGeo_
 
edm::ESGetToken< DTGeometry,
MuonGeometryRecord
dtGeomH_
 
int m_bx_window_
 
bool m_debug_
 
int m_max_quality_to_overwrite_t0_
 
double m_phi_window_
 
bool m_storeAllRPCHits_
 
RPCGeometry const * rpcGeo_
 
edm::ESGetToken< RPCGeometry,
MuonGeometryRecord
rpcGeomH_
 
double shift_back_
 

Static Private Attributes

static constexpr double distance_between_two_rpc_layers_ = 35
 
static constexpr double m_dt_phi_granularity_ = (65536. / 0.8)
 
static constexpr double m_dt_phiB_granularity_ = (2048. / 1.4)
 
static constexpr double R_ [2][2] = {{410.0, 444.8}, {492.7, 527.3}}
 

Detailed Description

Definition at line 56 of file RPCIntegrator.h.

Constructor & Destructor Documentation

RPCIntegrator::RPCIntegrator ( const edm::ParameterSet pset,
edm::ConsumesCollector iC 
)

Definition at line 13 of file RPCIntegrator.cc.

References edm::ConsumesCollector::esConsumes(), edm::ParameterSet::getUntrackedParameter(), and LogDebug.

13  {
14  m_debug_ = pset.getUntrackedParameter<bool>("debug");
15  if (m_debug_)
16  LogDebug("RPCIntegrator") << "RPCIntegrator constructor";
17  m_max_quality_to_overwrite_t0_ = pset.getUntrackedParameter<int>("max_quality_to_overwrite_t0");
18  m_bx_window_ = pset.getUntrackedParameter<int>("bx_window");
19  m_phi_window_ = pset.getUntrackedParameter<double>("phi_window");
20  m_storeAllRPCHits_ = pset.getUntrackedParameter<bool>("storeAllRPCHits");
21 
24 }
double m_phi_window_
Definition: RPCIntegrator.h:88
T getUntrackedParameter(std::string const &, T const &) const
int m_max_quality_to_overwrite_t0_
Definition: RPCIntegrator.h:86
bool m_storeAllRPCHits_
Definition: RPCIntegrator.h:89
edm::ESGetToken< DTGeometry, MuonGeometryRecord > dtGeomH_
Definition: RPCIntegrator.h:93
edm::ESGetToken< RPCGeometry, MuonGeometryRecord > rpcGeomH_
Definition: RPCIntegrator.h:94
#define LogDebug(id)
RPCIntegrator::~RPCIntegrator ( )

Definition at line 26 of file RPCIntegrator.cc.

References LogDebug.

26  {
27  if (m_debug_)
28  LogDebug("RPCIntegrator") << "RPCIntegrator destructor";
29 }
#define LogDebug(id)

Member Function Documentation

L1Phase2MuDTPhDigi RPCIntegrator::createL1Phase2MuDTPhDigi ( RPCDetId  rpcDetId,
int  rpc_bx,
double  rpc_time,
double  rpc_global_phi,
double  phiB,
int  rpc_flag 
)

Definition at line 189 of file RPCIntegrator.cc.

References RPCDetId::layer(), LogDebug, RPCDetId::ring(), RPCDetId::sector(), and RPCDetId::station().

190  {
191  if (m_debug_)
192  LogDebug("RPCIntegrator") << "Creating DT TP out of RPC recHits";
193  int rpc_wheel = rpcDetId.ring(); // In barrel, wheel is accessed via ring() method ([-2,+2])
194  int trigger_sector = rpcDetId.sector() - 1; // DT Trigger sector:[0,11] while RPC sector:[1,12]
195  int rpc_station = rpcDetId.station();
196  int rpc_layer = rpcDetId.layer();
197  int rpc_trigger_phi = phiInDTTPFormat(rpc_global_phi, rpcDetId.sector());
198  int rpc_trigger_phiB = (phiB == -10000) ? phiB : (int)round(phiB * m_dt_phiB_granularity_);
199  int rpc_quality = -1; // dummy for rpc
200  int rpc_index = 0; // dummy for rpc
201  return L1Phase2MuDTPhDigi(rpc_bx,
202  rpc_wheel,
203  trigger_sector,
204  rpc_station,
205  rpc_layer, //this would be the layer in the new dataformat
206  rpc_trigger_phi,
207  rpc_trigger_phiB,
208  rpc_quality,
209  rpc_index,
210  rpc_time,
211  -1, // no chi2 for RPC
212  rpc_flag);
213 }
static constexpr double m_dt_phiB_granularity_
Definition: RPCIntegrator.h:97
int ring() const
Definition: RPCDetId.h:59
int layer() const
Definition: RPCDetId.h:85
int sector() const
Sector id: the group of chambers at same phi (and increasing r)
Definition: RPCDetId.h:81
int phiInDTTPFormat(double rpc_global_phi, int rpcSector)
#define LogDebug(id)
int station() const
Definition: RPCDetId.h:78
void RPCIntegrator::finish ( )

Definition at line 44 of file RPCIntegrator.cc.

Referenced by progressbar.ProgressBar::__next__().

44 {}
bool RPCIntegrator::hasPosRF_rpc ( int  wh,
int  sec 
) const

Definition at line 257 of file RPCIntegrator.cc.

257 { return (wh > 0 || (wh == 0 && sec % 4 > 1)); }
void RPCIntegrator::initialise ( const edm::EventSetup iEventSetup,
double  shift_back_fromDT 
)

Definition at line 31 of file RPCIntegrator.cc.

References relativeConstraints::geom, edm::EventSetup::get(), edm::eventsetup::DependentRecordImplementation< RecordT, ListT >::get(), and LogDebug.

31  {
32  if (m_debug_)
33  LogDebug("RPCIntegrator") << "RPCIntegrator initialisation";
34 
35  if (m_debug_)
36  LogDebug("RPCIntegrator") << "Getting RPC geometry";
37 
38  const MuonGeometryRecord& geom = iEventSetup.get<MuonGeometryRecord>();
39  dtGeo_ = &geom.get(dtGeomH_);
40  rpcGeo_ = &geom.get(rpcGeomH_);
41  shift_back_ = shift_back_fromDT;
42 }
RPCGeometry const * rpcGeo_
Definition: RPCIntegrator.h:92
ProductT const & get(ESGetToken< ProductT, DepRecordT > const &iToken) const
double shift_back_
T get() const
Definition: EventSetup.h:88
edm::ESGetToken< DTGeometry, MuonGeometryRecord > dtGeomH_
Definition: RPCIntegrator.h:93
edm::ESGetToken< RPCGeometry, MuonGeometryRecord > rpcGeomH_
Definition: RPCIntegrator.h:94
#define LogDebug(id)
DTGeometry const * dtGeo_
Definition: RPCIntegrator.h:91
void RPCIntegrator::makeRPCOnlySegments ( )

Definition at line 76 of file RPCIntegrator.cc.

References funct::abs(), RPCRecHit::BunchX(), RPCDetId::layer(), SiStripPI::max, PV3DBase< T, PVType, FrameType >::phi(), RPCDetId::ring(), cmsdt::RPC_ASSOCIATE, RPCMetaprimitive::rpc_t0, RPCMetaprimitive::rpcFlag, RPCDetId::sector(), RPCDetId::station(), and FrontierCondition_GT_autoExpress_cfi::t0.

76  {
77  std::vector<L1Phase2MuDTPhDigi> rpc_only_segments;
78  for (auto& rpc_mp_it_layer1 : RPCMetaprimitives_) {
79  RPCDetId rpc_id_l1 = rpc_mp_it_layer1.rpc_id;
80  const RPCRecHit* rpc_cluster_l1 = rpc_mp_it_layer1.rpc_cluster;
81  GlobalPoint rpc_gp_l1 = rpc_mp_it_layer1.global_position;
82  if (rpc_id_l1.station() > 2 || rpc_id_l1.layer() != 1 ||
83  (rpc_mp_it_layer1.rpcFlag == RPC_ASSOCIATE && !m_storeAllRPCHits_))
84  continue;
85  // only one RPC layer in station three and four &&
86  // avoid duplicating pairs &&
87  // avoid building RPC only segment if DT segment was already there
88  int min_dPhi = std::numeric_limits<int>::max();
89  RPCMetaprimitive* bestMatch_rpc_mp_layer2 = nullptr;
90  for (auto& rpc_mp_it_layer2 : RPCMetaprimitives_) {
91  RPCDetId rpc_id_l2 = rpc_mp_it_layer2.rpc_id;
92  const RPCRecHit* rpc_cluster_l2 = rpc_mp_it_layer2.rpc_cluster;
93  GlobalPoint rpc_gp_l2 = rpc_mp_it_layer2.global_position;
94  if (rpc_id_l2.station() == rpc_id_l1.station() && rpc_id_l2.ring() == rpc_id_l1.ring() &&
95  rpc_id_l2.layer() != rpc_id_l1.layer() // ensure to have layer 1 --> layer 2
96  && rpc_id_l2.sector() == rpc_id_l1.sector() && rpc_cluster_l2->BunchX() == rpc_cluster_l1->BunchX() &&
97  (rpc_mp_it_layer2.rpcFlag != RPC_ASSOCIATE || m_storeAllRPCHits_)) {
98  // avoid building RPC only segment with a hit already matched to DT,
99  // except if one aske to store all RPC info
100  float tmp_dPhi = rpc_gp_l1.phi() - rpc_gp_l2.phi();
101  if (std::abs(tmp_dPhi) < std::abs(min_dPhi)) {
102  min_dPhi = tmp_dPhi;
103  bestMatch_rpc_mp_layer2 = &rpc_mp_it_layer2;
104  }
105  }
106  }
107  if (bestMatch_rpc_mp_layer2) {
108  rpc_mp_it_layer1.rpcFlag = 6;
109  // need a new flag (will be removed later) to differentiate
110  // between "has been matched to DT" and "Has been used in an RPC only segment"
111  bestMatch_rpc_mp_layer2->rpcFlag = 6;
112  double phiB = phiBending(&rpc_mp_it_layer1, &*bestMatch_rpc_mp_layer2);
113  // Arbitrarily choose the phi from layer 1
114  double global_phi = rpc_mp_it_layer1.global_position.phi();
115  double t0 = (rpc_mp_it_layer1.rpc_t0 + bestMatch_rpc_mp_layer2->rpc_t0) / 2;
116  // RPC only segment have rpcFlag==2
117  L1Phase2MuDTPhDigi rpc_only_segment =
118  createL1Phase2MuDTPhDigi(rpc_id_l1, rpc_mp_it_layer1.rpc_bx, t0, global_phi, phiB, 2);
119  rpc_only_segments.push_back(rpc_only_segment);
120  }
121  }
122  rpcRecHits_translated_.insert(rpcRecHits_translated_.end(), rpc_only_segments.begin(), rpc_only_segments.end());
123 }
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
L1Phase2MuDTPhDigi createL1Phase2MuDTPhDigi(RPCDetId rpcDetId, int rpc_bx, double rpc_time, double rpc_global_phi, double phiB, int rpc_flag)
std::vector< RPCMetaprimitive > RPCMetaprimitives_
Definition: RPCIntegrator.h:81
int ring() const
Definition: RPCDetId.h:59
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int layer() const
Definition: RPCDetId.h:85
std::vector< L1Phase2MuDTPhDigi > rpcRecHits_translated_
Definition: RPCIntegrator.h:80
int sector() const
Sector id: the group of chambers at same phi (and increasing r)
Definition: RPCDetId.h:81
double phiBending(RPCMetaprimitive *rpc_hit_1, RPCMetaprimitive *rpc_hit_2)
bool m_storeAllRPCHits_
Definition: RPCIntegrator.h:89
int BunchX() const
Definition: RPCRecHit.h:73
int station() const
Definition: RPCDetId.h:78
RPCMetaprimitive * RPCIntegrator::matchDTwithRPC ( cmsdt::metaPrimitive dt_metaprimitive)

Definition at line 153 of file RPCIntegrator.cc.

References funct::abs(), kinem::delta_phi(), DTChamberId, SiStripPI::max, cmsdt::metaPrimitive::phi, cmsdt::metaPrimitive::rawId, RPCDetId::ring(), cmsdt::RPC_ASSOCIATE, RPCMetaprimitive::rpcFlag, DTChamberId::sector(), RPCDetId::sector(), DTChamberId::station(), RPCDetId::station(), cmsdt::metaPrimitive::t0, and DTChamberId::wheel().

153  {
154  // metaprimitive dtChId is still in convention with [1 - 12]
155  // because at this stage the BX of metaprimitive is not yet computed...
156  // will also have to subtract 20*25 ns because of the recent change
157  int dt_bx = (int)round(dt_metaprimitive->t0 / 25.) - shift_back_;
158  DTChamberId dt_chId = DTChamberId(dt_metaprimitive->rawId);
159  int dt_sector = dt_chId.sector();
160  if (dt_sector == 13)
161  dt_sector = 4;
162  if (dt_sector == 14)
163  dt_sector = 10;
164  RPCMetaprimitive* bestMatch_rpcRecHit = nullptr;
165  float min_dPhi = std::numeric_limits<float>::max();
166  for (auto rpc_mp_it = RPCMetaprimitives_.begin(); rpc_mp_it != RPCMetaprimitives_.end(); rpc_mp_it++) {
167  RPCDetId rpc_det_id = rpc_mp_it->rpc_id;
168  if (rpc_det_id.ring() == dt_chId.wheel() // ring() in barrel RPC corresponds to the wheel
169  && rpc_det_id.station() == dt_chId.station() && rpc_det_id.sector() == dt_sector &&
170  std::abs(rpc_mp_it->rpc_bx - dt_bx) <= m_bx_window_) {
171  // Select the RPC hit closest in phi to the DT meta primitive
172 
173  // just a trick to apply the phi window cut on what could be accessed to fine tune it
174  int delta_phi =
175  (int)round((phi_DT_MP_conv(rpc_mp_it->global_position.phi(), rpc_det_id.sector()) - dt_metaprimitive->phi) *
177  if (std::abs(delta_phi) < min_dPhi && std::abs(delta_phi) < m_phi_window_) {
178  min_dPhi = std::abs(delta_phi);
179  bestMatch_rpcRecHit = &*rpc_mp_it;
180  }
181  }
182  }
183  if (bestMatch_rpcRecHit) {
184  bestMatch_rpcRecHit->rpcFlag = RPC_ASSOCIATE;
185  }
186  return bestMatch_rpcRecHit;
187 }
double m_phi_window_
Definition: RPCIntegrator.h:88
double phi_DT_MP_conv(double rpc_global_phi, int rpcSector)
std::vector< RPCMetaprimitive > RPCMetaprimitives_
Definition: RPCIntegrator.h:81
static constexpr double m_dt_phiB_granularity_
Definition: RPCIntegrator.h:97
int ring() const
Definition: RPCDetId.h:59
double shift_back_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double delta_phi(double ph11, double phi2)
Definition: AnglesUtil.h:84
int sector() const
Sector id: the group of chambers at same phi (and increasing r)
Definition: RPCDetId.h:81
int sector() const
Definition: DTChamberId.h:49
int station() const
Return the station number.
Definition: DTChamberId.h:42
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:39
int station() const
Definition: RPCDetId.h:78
void RPCIntegrator::matchWithDTAndUseRPCTime ( std::vector< cmsdt::metaPrimitive > &  dt_metaprimitives)

Definition at line 62 of file RPCIntegrator.cc.

References cmsdt::RPC_CONFIRM, RPCMetaprimitive::rpc_t0, and cmsdt::RPC_TIME.

62  {
63  for (auto dt_metaprimitive = dt_metaprimitives.begin(); dt_metaprimitive != dt_metaprimitives.end();
64  dt_metaprimitive++) {
65  RPCMetaprimitive* bestMatch_rpcRecHit = matchDTwithRPC(&*dt_metaprimitive);
66  if (bestMatch_rpcRecHit) {
67  (*dt_metaprimitive).rpcFlag = RPC_CONFIRM;
68  if ((*dt_metaprimitive).quality < m_max_quality_to_overwrite_t0_) {
69  (*dt_metaprimitive).t0 = bestMatch_rpcRecHit->rpc_t0 + 25 * shift_back_;
70  (*dt_metaprimitive).rpcFlag = RPC_TIME;
71  }
72  }
73  }
74 }
int m_max_quality_to_overwrite_t0_
Definition: RPCIntegrator.h:86
double shift_back_
RPCMetaprimitive * matchDTwithRPC(cmsdt::metaPrimitive *dt_metaprimitive)
double RPCIntegrator::phi_DT_MP_conv ( double  rpc_global_phi,
int  rpcSector 
)

Definition at line 235 of file RPCIntegrator.cc.

References edm::conversion(), and M_PI.

235  {
236  // Adaptation of https://github.com/cms-sw/cmssw/blob/master/L1Trigger/L1TTwinMux/src/RPCtoDTTranslator.cc#L349
237 
238  if (rpcSector == 1)
239  return rpc_global_phi;
240  else {
241  float conversion = 1 / 6.;
242  if (rpc_global_phi >= 0)
243  return rpc_global_phi - (rpcSector - 1) * M_PI * conversion;
244  else
245  return rpc_global_phi + (13 - rpcSector) * M_PI * conversion;
246  }
247 }
#define M_PI
void conversion(EventAux const &from, EventAuxiliary &to)
Definition: EventAux.cc:9
double RPCIntegrator::phiBending ( RPCMetaprimitive rpc_hit_1,
RPCMetaprimitive rpc_hit_2 
)

Definition at line 215 of file RPCIntegrator.cc.

References RPCMetaprimitive::global_position, PV3DBase< T, PVType, FrameType >::phi(), RPCDetId::ring(), RPCMetaprimitive::rpc_id, RPCDetId::sector(), slope, RPCDetId::station(), and PV3DBase< T, PVType, FrameType >::x().

215  {
216  DTChamberId DT_chamber(rpc_hit_1->rpc_id.ring(), rpc_hit_1->rpc_id.station(), rpc_hit_1->rpc_id.sector());
217  LocalPoint lp_rpc_hit_1_dtconv = dtGeo_->chamber(DT_chamber)->toLocal(rpc_hit_1->global_position);
218  LocalPoint lp_rpc_hit_2_dtconv = dtGeo_->chamber(DT_chamber)->toLocal(rpc_hit_2->global_position);
219  double slope = (lp_rpc_hit_1_dtconv.x() - lp_rpc_hit_2_dtconv.x()) / distance_between_two_rpc_layers_;
220  double average_x = (lp_rpc_hit_1_dtconv.x() + lp_rpc_hit_2_dtconv.x()) / 2;
221  GlobalPoint seg_middle_global =
222  dtGeo_->chamber(DT_chamber)->toGlobal(LocalPoint(average_x, 0., 0.)); // for station 1 and 2, z = 0
223  double seg_phi = phi_DT_MP_conv(seg_middle_global.phi(), rpc_hit_1->rpc_id.sector());
224  double psi = atan(slope);
225  double phiB = hasPosRF_rpc(rpc_hit_1->rpc_id.ring(), rpc_hit_1->rpc_id.sector()) ? psi - seg_phi : -psi - seg_phi;
226  return phiB;
227 }
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:30
const DTChamber * chamber(const DTChamberId &id) const
Return a DTChamber given its id.
Definition: DTGeometry.cc:90
static const double slope[3]
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
double phi_DT_MP_conv(double rpc_global_phi, int rpcSector)
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:58
GlobalPoint global_position
Definition: RPCIntegrator.h:38
std::map< std::string, int, std::less< std::string > > psi
int ring() const
Definition: RPCDetId.h:59
static constexpr double distance_between_two_rpc_layers_
bool hasPosRF_rpc(int wh, int sec) const
int sector() const
Sector id: the group of chambers at same phi (and increasing r)
Definition: RPCDetId.h:81
T x() const
Definition: PV3DBase.h:59
int station() const
Definition: RPCDetId.h:78
DTGeometry const * dtGeo_
Definition: RPCIntegrator.h:91
int RPCIntegrator::phiInDTTPFormat ( double  rpc_global_phi,
int  rpcSector 
)

Definition at line 229 of file RPCIntegrator.cc.

229  {
230  double rpc_localDT_phi;
231  rpc_localDT_phi = phi_DT_MP_conv(rpc_global_phi, rpcSector) * m_dt_phi_granularity_;
232  return (int)round(rpc_localDT_phi);
233 }
double phi_DT_MP_conv(double rpc_global_phi, int rpcSector)
static constexpr double m_dt_phi_granularity_
Definition: RPCIntegrator.h:96
void RPCIntegrator::prepareMetaPrimitives ( edm::Handle< RPCRecHitCollection rpcRecHits)

Definition at line 46 of file RPCIntegrator.cc.

References cmsdt::BX_SHIFT, cmsdt::LHC_CLK_FREQ, RPCDetId::region(), cmsdt::RPC_HIT, and RPCDetId.

46  {
47  RPCMetaprimitives_.clear();
48  rpcRecHits_translated_.clear();
49  for (const auto& rpcIt : *rpcRecHits) {
50  RPCDetId rpcDetId = (RPCDetId)(rpcIt).rpcId();
51  GlobalPoint global_position = RPCGlobalPosition(rpcDetId, rpcIt);
52  int rpc_region = rpcDetId.region();
53  if (rpc_region != 0)
54  continue; // Region = 0 Barrel
55 
56  // set everyone to rpc single hit (3) not matched to DT flag for now
57  // change last two elements if dt bx centered at zero again
58  RPCMetaprimitives_.emplace_back(
59  rpcDetId, &rpcIt, global_position, RPC_HIT, rpcIt.BunchX() + BX_SHIFT, rpcIt.time() + BX_SHIFT * LHC_CLK_FREQ);
60  }
61 }
std::vector< RPCMetaprimitive > RPCMetaprimitives_
Definition: RPCIntegrator.h:81
GlobalPoint RPCGlobalPosition(RPCDetId rpcId, const RPCRecHit &rpcIt) const
constexpr int BX_SHIFT
Definition: constants.h:257
std::vector< L1Phase2MuDTPhDigi > rpcRecHits_translated_
Definition: RPCIntegrator.h:80
constexpr int LHC_CLK_FREQ
Definition: constants.h:176
int region() const
Region id: 0 for Barrel, +/-1 For +/- Endcap.
Definition: RPCDetId.h:53
void RPCIntegrator::removeRPCHitsUsed ( )

Definition at line 136 of file RPCIntegrator.cc.

References LogDebug, and cmsdt::RPC_ASSOCIATE.

136  {
137  if (m_debug_)
138  LogDebug("RPCIntegrator") << "RPCIntegrator removeRPCHitsUsed method";
139  if (!m_storeAllRPCHits_) {
140  // Remove RPC hit attached to a DT or RPC segment if required by user
141  // (avoid having two TP's corresponding to the same physical hit)
142  auto rpcRecHit_translated_ = rpcRecHits_translated_.begin();
143  while (rpcRecHit_translated_ != rpcRecHits_translated_.end()) {
144  if (rpcRecHit_translated_->rpcFlag() == RPC_ASSOCIATE || rpcRecHit_translated_->rpcFlag() == 6) {
145  rpcRecHit_translated_ = rpcRecHits_translated_.erase(rpcRecHit_translated_);
146  } else {
147  ++rpcRecHit_translated_;
148  }
149  }
150  }
151 }
std::vector< L1Phase2MuDTPhDigi > rpcRecHits_translated_
Definition: RPCIntegrator.h:80
bool m_storeAllRPCHits_
Definition: RPCIntegrator.h:89
#define LogDebug(id)
GlobalPoint RPCIntegrator::RPCGlobalPosition ( RPCDetId  rpcId,
const RPCRecHit rpcIt 
) const

Definition at line 249 of file RPCIntegrator.cc.

References RPCRecHit::localPosition(), and RPCDetId.

249  {
250  RPCDetId rpcid = RPCDetId(rpcId);
251  const LocalPoint& rpc_lp = rpcIt.localPosition();
252  const GlobalPoint& rpc_gp = rpcGeo_->idToDet(rpcid)->surface().toGlobal(rpc_lp);
253 
254  return rpc_gp;
255 }
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:79
LocalPoint localPosition() const override
Return the 3-dimensional local position.
Definition: RPCRecHit.h:37
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
RPCGeometry const * rpcGeo_
Definition: RPCIntegrator.h:92
const GeomDet * idToDet(DetId) const override
Definition: RPCGeometry.cc:32
void RPCIntegrator::storeRPCSingleHits ( )

Definition at line 125 of file RPCIntegrator.cc.

References cmsdt::RPC_ASSOCIATE.

125  {
126  for (auto rpc_mp_it = RPCMetaprimitives_.begin(); rpc_mp_it != RPCMetaprimitives_.end(); rpc_mp_it++) {
127  RPCDetId rpcDetId = rpc_mp_it->rpc_id;
128  if (rpc_mp_it->rpcFlag == 6)
129  rpc_mp_it->rpcFlag = RPC_ASSOCIATE;
131  rpcDetId, rpc_mp_it->rpc_bx, rpc_mp_it->rpc_t0, rpc_mp_it->global_position.phi(), -10000, rpc_mp_it->rpcFlag);
132  rpcRecHits_translated_.push_back(rpc_out);
133  }
134 }
L1Phase2MuDTPhDigi createL1Phase2MuDTPhDigi(RPCDetId rpcDetId, int rpc_bx, double rpc_time, double rpc_global_phi, double phiB, int rpc_flag)
std::vector< RPCMetaprimitive > RPCMetaprimitives_
Definition: RPCIntegrator.h:81
std::vector< L1Phase2MuDTPhDigi > rpcRecHits_translated_
Definition: RPCIntegrator.h:80

Member Data Documentation

constexpr double RPCIntegrator::distance_between_two_rpc_layers_ = 35
staticprivate

Definition at line 101 of file RPCIntegrator.h.

DTGeometry const* RPCIntegrator::dtGeo_
private

Definition at line 91 of file RPCIntegrator.h.

edm::ESGetToken<DTGeometry, MuonGeometryRecord> RPCIntegrator::dtGeomH_
private

Definition at line 93 of file RPCIntegrator.h.

int RPCIntegrator::m_bx_window_
private

Definition at line 87 of file RPCIntegrator.h.

bool RPCIntegrator::m_debug_
private

Definition at line 85 of file RPCIntegrator.h.

constexpr double RPCIntegrator::m_dt_phi_granularity_ = (65536. / 0.8)
staticprivate

Definition at line 96 of file RPCIntegrator.h.

constexpr double RPCIntegrator::m_dt_phiB_granularity_ = (2048. / 1.4)
staticprivate

Definition at line 97 of file RPCIntegrator.h.

int RPCIntegrator::m_max_quality_to_overwrite_t0_
private

Definition at line 86 of file RPCIntegrator.h.

double RPCIntegrator::m_phi_window_
private

Definition at line 88 of file RPCIntegrator.h.

bool RPCIntegrator::m_storeAllRPCHits_
private

Definition at line 89 of file RPCIntegrator.h.

constexpr double RPCIntegrator::R_[2][2] = {{410.0, 444.8}, {492.7, 527.3}}
staticprivate

Definition at line 100 of file RPCIntegrator.h.

RPCGeometry const* RPCIntegrator::rpcGeo_
private

Definition at line 92 of file RPCIntegrator.h.

edm::ESGetToken<RPCGeometry, MuonGeometryRecord> RPCIntegrator::rpcGeomH_
private

Definition at line 94 of file RPCIntegrator.h.

std::vector<RPCMetaprimitive> RPCIntegrator::RPCMetaprimitives_

Definition at line 81 of file RPCIntegrator.h.

std::vector<L1Phase2MuDTPhDigi> RPCIntegrator::rpcRecHits_translated_

Definition at line 80 of file RPCIntegrator.h.

double RPCIntegrator::shift_back_
private

Definition at line 103 of file RPCIntegrator.h.