CMS 3D CMS Logo

List of all members | Public Member Functions
CPPFMaskReClusterizer Class Reference

#include <CPPFMaskReClusterizer.h>

Public Member Functions

 CPPFMaskReClusterizer ()
 
CPPFClusterContainer doAction (const RPCDetId &id, CPPFClusterContainer &initClusters, const CPPFRollMask &mask) const
 
bool get (const CPPFRollMask &mask, int strip) const
 
 ~CPPFMaskReClusterizer ()
 

Detailed Description

\Class CPPFMaskReClusterizer

Author
R. Hadjiiska – INRNE-BAS, Sofia

Definition at line 12 of file CPPFMaskReClusterizer.h.

Constructor & Destructor Documentation

◆ CPPFMaskReClusterizer()

CPPFMaskReClusterizer::CPPFMaskReClusterizer ( )
inline

Definition at line 14 of file CPPFMaskReClusterizer.h.

14 {};

◆ ~CPPFMaskReClusterizer()

CPPFMaskReClusterizer::~CPPFMaskReClusterizer ( )
inline

Definition at line 15 of file CPPFMaskReClusterizer.h.

15 {};

Member Function Documentation

◆ doAction()

CPPFClusterContainer CPPFMaskReClusterizer::doAction ( const RPCDetId id,
CPPFClusterContainer initClusters,
const CPPFRollMask mask 
) const

\Class CPPFMaskReClusterizer

Author
R. Hadjiiska – INRNE-BAS, Sofia

Definition at line 9 of file CPPFMaskReClusterizer.cc.

11  {
12  CPPFClusterContainer finClusters;
13  if (initClusters.empty())
14  return finClusters;
15 
16  CPPFCluster prev = *initClusters.begin();
17  for (auto cl = std::next(initClusters.begin()); cl != initClusters.end(); ++cl) {
18  // Merge this cluster if it is adjacent by 1 masked strip
19  // Note that the CPPFClusterContainer collection is sorted in DECREASING ORDER of strip #
20  // So the prev. cluster is placed after the current cluster (check the < operator of CPPFCluster carefully)
21  if ((prev.firstStrip() - cl->lastStrip()) == 2 and this->get(mask, cl->lastStrip() + 1) and prev.bx() == cl->bx()) {
22  CPPFCluster merged(cl->firstStrip(), prev.lastStrip(), cl->bx());
23  prev = merged;
24  } else {
25  finClusters.insert(prev);
26  prev = *cl;
27  }
28  }
29 
30  // Finalize by putting the last cluster to the collection
31  finClusters.insert(prev);
32 
33  return finClusters;
34 }

References CPPFCluster::bx(), GetRecoTauVFromDQM_MC_cff::cl, CPPFCluster::firstStrip(), get(), CPPFCluster::lastStrip(), and GetRecoTauVFromDQM_MC_cff::next.

Referenced by RecHitProcessor::process(), and RecHitProcessor::processLook().

◆ get()

bool CPPFMaskReClusterizer::get ( const CPPFRollMask mask,
int  strip 
) const
CPPFClusterContainer
std::set< CPPFCluster > CPPFClusterContainer
Definition: CPPFClusterContainer.h:4
CPPFCluster::firstStrip
int firstStrip() const
Definition: CPPFCluster.cc:16
digitizers_cfi.strip
strip
Definition: digitizers_cfi.py:19
GetRecoTauVFromDQM_MC_cff.cl
cl
Definition: GetRecoTauVFromDQM_MC_cff.py:38
CPPFCluster::bx
int bx() const
Definition: CPPFCluster.cc:19
CPPFMaskReClusterizer::get
bool get(const CPPFRollMask &mask, int strip) const
Definition: CPPFMaskReClusterizer.cc:36
CPPFCluster::lastStrip
int lastStrip() const
Definition: CPPFCluster.cc:17
GetRecoTauVFromDQM_MC_cff.next
next
Definition: GetRecoTauVFromDQM_MC_cff.py:31
CPPFCluster
Definition: CPPFCluster.h:4