22 #include <Math/Functions.h>
23 #include <Math/SVector.h>
24 #include <Math/SMatrix.h>
58 LogDebug(
"GE0SegAlgoRU") <<
myName <<
" has algorithm cuts set to: \n"
59 <<
"--------------------------------------------------------------------\n"
61 <<
"doCollisions = " << doCollisions <<
"\n"
75 auto& superchamber = ensemble.first;
76 for (
const auto& rechit : rechits) {
80 hitAndPositions.emplace_back(&(*rechit), nLoc, glb, hitAndPositions.size());
83 LogDebug(
"GE0Segment|GE0") <<
"found " << hitAndPositions.size() <<
" rechits in superchamber " << superchamber->id();
85 float z1 = superchamber->chamber(1)->position().z();
86 float zlast = superchamber->chamber(superchamber->nChambers())->
position().z();
89 return h1.
layer < h2.layer;
93 return h1.
layer > h2.layer;
95 return run(ensemble.first, hitAndPositions);
99 #ifdef EDM_ML_DEBUG // have lines below only compiled when in debug mode
101 edm::LogVerbatim(
"GE0SegAlgoRU") <<
"[GEMSegmentAlgorithm::run] build segments in chamber " << chId
102 <<
" which contains " << rechits.size() <<
" rechits";
103 for (
const auto&
h : rechits) {
104 auto ge0id =
h.rh->gemId();
106 edm::LogVerbatim(
"GE0SegAlgoRU") <<
"[RecHit :: Loc x = " << std::showpos << std::setw(9) << rhLP.x()
107 <<
" Glb y = " << std::showpos << std::setw(9) << rhLP.y()
108 <<
" Time = " << std::showpos <<
h.rh->BunchX() <<
" -- " << ge0id.rawId() <<
" = "
109 << ge0id <<
" ]" << std::endl;
114 return std::vector<GEMSegment>();
120 for (
const auto& rechit : rechits) {
121 recHits_per_layer[rechit.layer - 1]++;
140 std::vector<GEMSegment> segments;
146 auto doDisplaced = [&]() {
155 auto printSegments = [&] {
156 #ifdef EDM_ML_DEBUG // have lines below only compiled when in debug mode
157 for (
const auto& seg : segments) {
159 const auto& rechits = seg.specificRecHits();
160 edm::LogVerbatim(
"GE0SegAlgoRU") <<
"[GE0SegAlgoRU] segment in chamber " << chId <<
" which contains "
161 << rechits.size() <<
" rechits and with specs: \n"
163 for (
const auto& rh : rechits) {
164 auto ge0id = rh.gemId();
165 edm::LogVerbatim(
"GE0SegAlgoRU") <<
"[RecHit :: Loc x = " << std::showpos << std::setw(9)
166 << rh.localPosition().x() <<
" Loc y = " << std::showpos << std::setw(9)
167 << rh.localPosition().y() <<
" Time = " << std::showpos << rh.BunchX()
168 <<
" -- " << ge0id.rawId() <<
" = " << ge0id <<
" ]";
201 const unsigned int n_seg_min,
203 const std::vector<unsigned int>& recHits_per_layer,
205 std::vector<GEMSegment>& segments)
const {
206 auto ib = rechits.begin();
207 auto ie = rechits.end();
208 std::vector<std::pair<float, HitAndPositionPtrContainer> > proto_segments;
210 for (
auto i1 =
ib; i1 != ie; ++i1) {
211 const auto& h1 = *i1;
220 for (
auto i2 = ie - 1; i2 != i1; --i2) {
221 const auto& h2 = *i2;
230 if ((
std::abs(
int(h2.layer) -
int(h1.layer)) + 1) <
int(n_seg_min))
239 std::unique_ptr<MuonSegFit> current_fit;
240 current_fit =
addHit(current_proto_segment, h1);
241 current_fit =
addHit(current_proto_segment, h2);
247 current_proto_segment,
252 if (current_proto_segment.size() > n_seg_min)
255 edm::LogVerbatim(
"GE0SegAlgoRU") <<
"[GE0SegAlgoRU::lookForSegments] # of hits in segment "
256 << current_proto_segment.size() <<
" min # " << n_seg_min <<
" => "
257 << (current_proto_segment.size() >= n_seg_min) <<
" chi2/ndof "
258 << current_fit->chi2() / current_fit->ndof() <<
" => "
259 << (current_fit->chi2() / current_fit->ndof() < params.
maxChi2GoodSeg)
262 if (current_proto_segment.size() < n_seg_min)
264 const float current_metric = current_fit->chi2() / current_fit->ndof();
271 for (
const auto* rh : current_proto_segment) {
277 if (nNonCentral >= nCentral)
281 proto_segments.emplace_back(current_metric, current_proto_segment);
288 std::vector<GEMSegment>& segments,
290 std::sort(proto_segments.begin(),
291 proto_segments.end(),
292 [](
const std::pair<float, HitAndPositionPtrContainer>&
a,
293 const std::pair<float, HitAndPositionPtrContainer>&
b) {
return a.first <
b.first; });
296 std::vector<unsigned int> usedHits;
297 for (
auto& container : proto_segments) {
301 bool alreadyFilled =
false;
302 for (
const auto&
h : currentProtoSegment) {
303 for (
unsigned int iOH = 0; iOH < usedHits.size(); ++iOH) {
304 if (usedHits[iOH] !=
h->idx)
306 alreadyFilled =
true;
312 for (
const auto*
h : currentProtoSegment) {
313 usedHits.push_back(
h->idx);
317 std::unique_ptr<MuonSegFit> current_fit =
makeFit(currentProtoSegment);
322 float averageBX = 0.;
323 for (
const auto*
h : currentProtoSegment) {
324 averageBX +=
h->rh->BunchX();
326 averageBX /= int(currentProtoSegment.size());
328 std::sort(currentProtoSegment.begin(),
329 currentProtoSegment.end(),
332 std::vector<const GEMRecHit*> bareRHs;
333 bareRHs.reserve(currentProtoSegment.size());
334 for (
const auto* rh : currentProtoSegment)
335 bareRHs.push_back(rh->rh);
338 current_fit->intercept(),
339 current_fit->localdir(),
340 current_fit->covarianceMatrix(),
344 segments.push_back(
temp);
351 std::unique_ptr<MuonSegFit>& current_fit,
354 HitAndPositionContainer::const_iterator i1,
355 HitAndPositionContainer::const_iterator i2)
const {
368 for (
auto iH = i1 + 1; iH != i2; ++iH) {
369 if (iH->layer == i1->layer)
371 if (iH->layer == i2->layer)
385 const bool beamConst,
389 edm::LogVerbatim(
"GE0SegAlgoRU") <<
"[GE0SegAlgoRU::areHitsCloseInEta] gp1 = " << h1 <<
" in eta part = " << h1.
eta()
390 <<
" and gp2 = " << h2 <<
" in eta part = " << h2.
eta() <<
" ==> dEta = " << diff
391 <<
" ==> return " << (diff < 0.1) << std::endl;
397 const unsigned int nLayDisp,
401 edm::LogVerbatim(
"GE0SegAlgoRU") <<
"[GE0SegAlgoRU::areHitsCloseInGlobalPhi] gp1 = " << h1 <<
" and gp2 = " << h2
402 <<
" ==> dPhi = " << dphi12 <<
" ==> return "
409 const std::unique_ptr<MuonSegFit>&
fit,
413 const float avgETA = (proto_segment[1]->gp.eta() + proto_segment[0]->gp.eta()) / 2.;
425 float x = fit->xfit(z);
426 float y = fit->yfit(z);
432 proto_segment.push_back(&aHit);
441 for (
const auto& rh : proto_segment) {
445 muonRecHits.push_back(trkRecHit);
447 auto currentFit = std::make_unique<MuonSegFit>(muonRecHits);
454 std::unique_ptr<MuonSegFit>&
fit,
455 const unsigned int n_seg_min)
const {
456 while (proto_segment.size() > n_seg_min && fit->chi2() / fit->ndof() >
maxChi2) {
458 HitAndPositionPtrContainer::iterator worstHit;
459 for (
auto it = proto_segment.begin(); it != proto_segment.end(); ++it) {
466 edm::LogVerbatim(
"GE0SegAlgoRU") <<
"[GE0SegAlgoRU::pruneBadHits] pruning one hit-> layer: " << (*worstHit)->layer
467 <<
" eta: " << (*worstHit)->gp.eta() <<
" phi: " << (*worstHit)->gp.phi()
468 <<
" old chi2/dof: " << fit->chi2() / fit->ndof() << std::endl;
469 proto_segment.erase(worstHit);
477 const float du = fit->xdev(lp.x(), lp.z());
478 const float dv = fit->ydev(lp.y(), lp.z());
480 ROOT::Math::SMatrix<double, 2, 2, ROOT::Math::MatRepSym<double, 2> > IC;
487 return du * du * IC(0, 0) + 2. * du * dv * IC(0, 1) + dv * dv * IC(1, 1);
491 for (
const auto*
h : proto_segment)
492 if (
h->layer == layer)
503 for (
auto it = new_proto_segment.begin(); it != new_proto_segment.end();) {
504 if ((*it)->layer == new_hit.
layer) {
506 it = new_proto_segment.erase(it);
511 if (old_hit ==
nullptr)
513 auto new_fit =
addHit(new_proto_segment, new_hit);
517 if (old_hit->
lp == new_hit.
lp) {
519 for (
const auto*
h : current_proto_segment)
521 avgbx +=
h->rh->BunchX();
522 avgbx /= float(current_proto_segment.size() - 1);
526 else if (new_fit->chi2() < current_fit->chi2())
530 current_proto_segment = new_proto_segment;
536 std::unique_ptr<MuonSegFit>& current_fit,
540 auto new_fit =
addHit(new_proto_segment, new_hit);
541 if (new_fit->chi2() / new_fit->ndof() <
maxChi2) {
542 current_proto_segment = new_proto_segment;
std::unique_ptr< MuonSegFit > addHit(HitAndPositionPtrContainer &proto_segment, const HitAndPosition &aHit) const
std::shared_ptr< TrackingRecHit > MuonRecHitPtr
Log< level::Info, true > LogVerbatim
std::vector< bool > BoolContainer
Point3DBase< Scalar, LocalTag > LocalPoint
float computeDeltaPhi(const LocalPoint &position, const LocalVector &direction) const
GE0SegAlgoRU(const edm::ParameterSet &ps)
Constructor.
SegmentParameters wideParameters
LocalPoint localPosition() const override
Return the 3-dimensional local position.
void increaseProtoSegment(const float maxChi2, std::unique_ptr< MuonSegFit > ¤t_fit, HitAndPositionPtrContainer ¤t_proto_segment, const HitAndPosition &new_hit) const
float getHitSegChi2(const std::unique_ptr< MuonSegFit > &fit, const GEMRecHit &hit) const
bool areHitsCloseInEta(const float maxETA, const bool beamConst, const GlobalPoint &h1, const GlobalPoint &h2) const
unsigned int maxNumberOfHitsPerLayer
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
std::vector< const HitAndPosition * > HitAndPositionPtrContainer
Geom::Phi< T > phi() const
std::pair< const GEMSuperChamber *, std::map< uint32_t, const GEMEtaPartition * > > GEMEnsemble
const GEMSuperChamber * theChamber
unsigned int maxNumberOfHits
void tryAddingHitsToSegment(const float maxETA, const float maxPhi, const float maxChi2, std::unique_ptr< MuonSegFit > ¤t_fit, HitAndPositionPtrContainer &proto_segment, const BoolContainer &used, HitAndPositionContainer::const_iterator i1, HitAndPositionContainer::const_iterator i2) const
constexpr std::array< uint8_t, layerIndexSize > layer
void compareProtoSegment(std::unique_ptr< MuonSegFit > ¤t_fit, HitAndPositionPtrContainer ¤t_proto_segment, const HitAndPosition &new_hit) const
void pruneBadHits(const float maxChi2, HitAndPositionPtrContainer &proto_segment, std::unique_ptr< MuonSegFit > &fit, const unsigned int n_seg_min) const
bool hasHitOnLayer(const HitAndPositionPtrContainer &proto_segment, const unsigned int layer) const
void lookForSegments(const SegmentParameters ¶ms, const unsigned int n_seg_min, const HitAndPositionContainer &rechits, const std::vector< unsigned int > &recHits_per_layer, BoolContainer &used, std::vector< GEMSegment > &segments) const
GEMRecHit * clone() const override
std::vector< MuonRecHitPtr > MuonRecHitContainer
LocalError localPositionError() const override
Return the 3-dimensional error on the local position.
Abs< T >::type abs(const T &t)
std::vector< GEMSegment > run(const GEMSuperChamber *chamber, const HitAndPositionContainer &rechits)
SegmentParameters displacedParameters
GlobalPoint globalAtZ(const std::unique_ptr< MuonSegFit > &fit, float z) const
unsigned int minNumberOfHits
int nChambers() const
Return numbers of chambers.
GEMDetId id() const
Return the GEMDetId of this super chamber.
T getParameter(std::string const &) const
void addUniqueSegments(SegmentByMetricContainer &proto_segments, std::vector< GEMSegment > &segments, BoolContainer &used) const
bool areHitsCloseInGlobalPhi(const float maxPHI, const unsigned int nLayDisp, const GlobalPoint &h1, const GlobalPoint &h2) const
std::vector< std::pair< float, HitAndPositionPtrContainer > > SegmentByMetricContainer
bool isHitNearSegment(const float maxETA, const float maxPHI, const std::unique_ptr< MuonSegFit > &fit, const HitAndPositionPtrContainer &proto_segment, const HitAndPosition &h) const
std::vector< HitAndPosition > HitAndPositionContainer
static int position[264][3]
SegmentParameters stdParameters
void setPosition(LocalPoint pos)
Set local position.
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
std::unique_ptr< MuonSegFit > makeFit(const HitAndPositionPtrContainer &proto_segment) const