CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Private Member Functions | Private Attributes | Friends
PFBlockAlgo Class Reference

Particle Flow Algorithm. More...

#include <PFBlockAlgo.h>

Public Types

typedef std::vector< std::unique_ptr< reco::PFBlockElement > > ElementList
 
typedef std::array< std::pair< unsigned int, unsigned int >, reco::PFBlockElement::kNBETypesElementRanges
 

Public Member Functions

void buildElements (const edm::Event &)
 
reco::PFBlockCollection findBlocks ()
 build blocks More...
 
 PFBlockAlgo ()
 
void setDebug (bool debug)
 sets debug printout flag More...
 
void setImporters (const std::vector< edm::ParameterSet > &, edm::ConsumesCollector &)
 
void setLinkers (const std::vector< edm::ParameterSet > &)
 
void updateEventSetup (const edm::EventSetup &)
 
 ~PFBlockAlgo ()
 

Private Member Functions

void link (const reco::PFBlockElement *el1, const reco::PFBlockElement *el2, double &dist) const
 check whether 2 elements are linked. Returns distance More...
 
void packLinks (reco::PFBlock &block, const std::unordered_map< std::pair< unsigned int, unsigned int >, double > &links) const
 

Private Attributes

bool debug_
 if true, debug printouts activated More...
 
ElementList elements_
 
const std::unordered_map< std::string, reco::PFBlockElement::TypeelementTypes_
 
std::vector< std::unique_ptr< BlockElementImporterBase > > importers_
 
std::vector< std::unique_ptr< KDTreeLinkerBase > > kdtrees_
 
std::vector< std::unique_ptr< BlockElementLinkerBase > > linkTests_
 
unsigned int linkTestSquare_ [reco::PFBlockElement::kNBETypes][reco::PFBlockElement::kNBETypes]
 
ElementRanges ranges_
 
bool useHO_
 

Friends

std::ostream & operator<< (std::ostream &, const PFBlockAlgo &)
 

Detailed Description

Particle Flow Algorithm.

Author
Colin Bernet (rewrite/refactor by L. Gray)
Date
January 2006 (April 2014)

Definition at line 34 of file PFBlockAlgo.h.

Member Typedef Documentation

◆ ElementList

typedef std::vector<std::unique_ptr<reco::PFBlockElement> > PFBlockAlgo::ElementList

Definition at line 37 of file PFBlockAlgo.h.

◆ ElementRanges

typedef std::array<std::pair<unsigned int, unsigned int>, reco::PFBlockElement::kNBETypes> PFBlockAlgo::ElementRanges

Definition at line 39 of file PFBlockAlgo.h.

Constructor & Destructor Documentation

◆ PFBlockAlgo()

PFBlockAlgo::PFBlockAlgo ( )

Definition at line 69 of file PFBlockAlgo.cc.

References ECAL, HCAL, DigiToRawDM_cff::HO, and INIT_ENTRY.

70  : debug_(false),
71  elementTypes_({INIT_ENTRY(PFBlockElement::TRACK),
72  INIT_ENTRY(PFBlockElement::PS1),
73  INIT_ENTRY(PFBlockElement::PS2),
76  INIT_ENTRY(PFBlockElement::GSF),
77  INIT_ENTRY(PFBlockElement::BREM),
78  INIT_ENTRY(PFBlockElement::HFEM),
79  INIT_ENTRY(PFBlockElement::HFHAD),
80  INIT_ENTRY(PFBlockElement::SC),
82  INIT_ENTRY(PFBlockElement::HGCAL)}) {}
#define INIT_ENTRY(name)
Definition: PFBlockAlgo.cc:17
const std::unordered_map< std::string, reco::PFBlockElement::Type > elementTypes_
Definition: PFBlockAlgo.h:82
bool debug_
if true, debug printouts activated
Definition: PFBlockAlgo.h:75

◆ ~PFBlockAlgo()

PFBlockAlgo::~PFBlockAlgo ( )

Definition at line 132 of file PFBlockAlgo.cc.

References gather_cfg::cout, debug_, and elements_.

132  {
133 #ifdef PFLOW_DEBUG
134  if (debug_)
135  cout << "~PFBlockAlgo - number of remaining elements: " << elements_.size() << endl;
136 #endif
137 }
ElementList elements_
Definition: PFBlockAlgo.h:71
bool debug_
if true, debug printouts activated
Definition: PFBlockAlgo.h:75

Member Function Documentation

◆ buildElements()

void PFBlockAlgo::buildElements ( const edm::Event evt)

Definition at line 284 of file PFBlockAlgo.cc.

References a, b, elements_, mps_fire::i, importers_, ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::it, kdtrees_, ranges_, and jetUpdater_cfi::sort.

Referenced by PFBlockProducer::produce().

284  {
285  // import block elements as defined in python configuration
286  ranges_.fill(std::make_pair(0, 0));
287  elements_.clear();
288  for (const auto& importer : importers_) {
289  importer->importToBlock(evt, elements_);
290  }
291 
292  std::sort(elements_.begin(), elements_.end(), [](const auto& a, const auto& b) { return a->type() < b->type(); });
293 
294  // list is now partitioned, so mark the boundaries so we can efficiently skip chunks
295  unsigned current_type = (!elements_.empty() ? elements_[0]->type() : 0);
296  unsigned last_type = (!elements_.empty() ? elements_.back()->type() : 0);
297  ranges_[current_type].first = 0;
298  ranges_[last_type].second = elements_.size() - 1;
299  for (size_t i = 0; i < elements_.size(); ++i) {
300  const auto the_type = elements_[i]->type();
301  if (the_type != current_type) {
302  ranges_[the_type].first = i;
303  ranges_[current_type].second = i - 1;
304  current_type = the_type;
305  }
306  }
307  // -------------- Loop over block elements ---------------------
308 
309  // Here we provide to all KDTree linkers the collections to link.
310  // Glowinski & Gouzevitch
311 
312  for (ElementList::iterator it = elements_.begin(); it != elements_.end(); ++it) {
313  for (const auto& kdtree : kdtrees_) {
314  if ((*it)->type() == kdtree->targetType()) {
315  kdtree->insertTargetElt(it->get());
316  }
317  if ((*it)->type() == kdtree->fieldType()) {
318  kdtree->insertFieldClusterElt(it->get());
319  }
320  }
321  }
322  //std::cout << "(new) imported: " << elements_.size() << " elements!" << std::endl;
323 }
ElementList elements_
Definition: PFBlockAlgo.h:71
std::vector< std::unique_ptr< KDTreeLinkerBase > > kdtrees_
Definition: PFBlockAlgo.h:86
ElementRanges ranges_
Definition: PFBlockAlgo.h:72
std::vector< std::unique_ptr< BlockElementImporterBase > > importers_
Definition: PFBlockAlgo.h:80
double b
Definition: hdecay.h:120
double a
Definition: hdecay.h:121

◆ findBlocks()

reco::PFBlockCollection PFBlockAlgo::findBlocks ( )

build blocks

Definition at line 139 of file PFBlockAlgo.cc.

References OfflinePixel3DPrimaryVertices_cfi::block_size, gather_cfg::blocks, elements_, mps_fire::i, dqmiolumiharvest::j, kdtrees_, submitPVResolutionJobs::key, relativeConstraints::keys, electronStore::links, linkTests_, linkTestSquare_, pfDeepBoostedJetPreprocessParams_cfi::lower_bound, SiStripPI::min, LaserDQM_cfg::p1, SiStripOfflineCRack_cfg::p2, packLinks(), FastTimerService_cff::range, and ranges_.

Referenced by PFBlockProducer::produce().

139  {
140  // Glowinski & Gouzevitch
141  for (const auto& kdtree : kdtrees_) {
142  kdtree->process();
143  }
144  // !Glowinski & Gouzevitch
146  // the blocks have not been passed to the event, and need to be cleared
147  blocks.reserve(elements_.size());
148 
149  QuickUnion qu(elements_.size());
150  const auto elem_size = elements_.size();
151  for (unsigned i = 0; i < elem_size; ++i) {
152  for (unsigned j = i + 1; j < elem_size; ++j) {
153  if (qu.connected(i, j))
154  continue;
155  if (!linkTests_[linkTestSquare_[elements_[i]->type()][elements_[j]->type()]]) {
156  j = ranges_[elements_[j]->type()].second;
157  continue;
158  }
159  auto p1(elements_[i].get()), p2(elements_[j].get());
160  const PFBlockElement::Type type1 = p1->type();
161  const PFBlockElement::Type type2 = p2->type();
162  const unsigned index = linkTestSquare_[type1][type2];
163  if (linkTests_[index]->linkPrefilter(p1, p2)) {
164  const double dist = linkTests_[index]->testLink(p1, p2);
165  // compute linking info if it is possible
166  if (dist > -0.5) {
167  qu.unite(i, j);
168  }
169  }
170  }
171  }
172 
173  std::unordered_multimap<unsigned, unsigned> blocksmap(elements_.size());
174  std::vector<unsigned> keys;
175  keys.reserve(elements_.size());
176  for (unsigned i = 0; i < elements_.size(); ++i) {
177  unsigned key = i;
178  while (key != qu.find(key))
179  key = qu.find(key); // make sure we always find the root node...
180  auto pos = std::lower_bound(keys.begin(), keys.end(), key);
181  if (pos == keys.end() || *pos != key) {
182  keys.insert(pos, key);
183  }
184  blocksmap.emplace(key, i);
185  }
186 
187  for (auto key : keys) {
188  blocks.push_back(reco::PFBlock());
189  auto range = blocksmap.equal_range(key);
190  auto& the_block = blocks.back();
191  ElementList::value_type::pointer p1(elements_[range.first->second].get());
192  the_block.addElement(p1);
193  const unsigned block_size = blocksmap.count(key) + 1;
194  //reserve up to 1M or 8MB; pay rehash cost for more
195  std::unordered_map<std::pair<unsigned int, unsigned int>, double> links(min(1000000u, block_size * block_size));
196  auto itr = range.first;
197  ++itr;
198  for (; itr != range.second; ++itr) {
199  ElementList::value_type::pointer p2(elements_[itr->second].get());
200  const PFBlockElement::Type type1 = p1->type();
201  const PFBlockElement::Type type2 = p2->type();
202  the_block.addElement(p2);
203  const unsigned index = linkTestSquare_[type1][type2];
204  if (nullptr != linkTests_[index]) {
205  const double dist = linkTests_[index]->testLink(p1, p2);
206  links.emplace(std::make_pair(p1->index(), p2->index()), dist);
207  }
208  }
209  packLinks(the_block, links);
210  }
211 
212  elements_.clear();
213 
214  return blocks;
215 }
ElementList elements_
Definition: PFBlockAlgo.h:71
std::vector< std::unique_ptr< KDTreeLinkerBase > > kdtrees_
Definition: PFBlockAlgo.h:86
key
prepare the HTCondor submission files and eventually submit them
void packLinks(reco::PFBlock &block, const std::unordered_map< std::pair< unsigned int, unsigned int >, double > &links) const
Definition: PFBlockAlgo.cc:217
std::vector< PFBlock > PFBlockCollection
collection of PFBlock objects
Definition: PFBlockFwd.h:10
ElementRanges ranges_
Definition: PFBlockAlgo.h:72
std::vector< std::unique_ptr< BlockElementLinkerBase > > linkTests_
Definition: PFBlockAlgo.h:83
unsigned int linkTestSquare_[reco::PFBlockElement::kNBETypes][reco::PFBlockElement::kNBETypes]
Definition: PFBlockAlgo.h:84
Block of elements.
Definition: PFBlock.h:26

◆ link()

void PFBlockAlgo::link ( const reco::PFBlockElement el1,
const reco::PFBlockElement el2,
double &  dist 
) const
inlineprivate

check whether 2 elements are linked. Returns distance

Definition at line 261 of file PFBlockAlgo.cc.

References ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), gather_cfg::cout, debug_, reco::PFBlockElement::kNBETypes, linkTests_, SiStripPI::max, SiStripPI::min, and reco::PFBlockElement::type().

Referenced by packLinks().

261  {
262  constexpr unsigned rowsize = reco::PFBlockElement::kNBETypes;
263  dist = -1.0;
264  const PFBlockElement::Type type1 = el1->type();
265  const PFBlockElement::Type type2 = el2->type();
266  const unsigned index = rowsize * std::max(type1, type2) + std::min(type1, type2);
267  if (debug_) {
268  std::cout << " PFBlockAlgo links type1 " << type1 << " type2 " << type2 << std::endl;
269  }
270 
271  // index is always checked in the preFilter above, no need to check here
272  dist = linkTests_[index]->testLink(el1, el2);
273 }
std::vector< std::unique_ptr< BlockElementLinkerBase > > linkTests_
Definition: PFBlockAlgo.h:83
bool debug_
if true, debug printouts activated
Definition: PFBlockAlgo.h:75

◆ packLinks()

void PFBlockAlgo::packLinks ( reco::PFBlock block,
const std::unordered_map< std::pair< unsigned int, unsigned int >, double > &  links 
) const
private

compute missing links in the blocks (the recursive procedure does not build all links)

Definition at line 217 of file PFBlockAlgo.cc.

References groupFilesInBlocks::block, ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), funct::false, testProducerWithPsetDescEmpty_cfi::i1, testProducerWithPsetDescEmpty_cfi::i2, reco::PFBlockElement::kNBETypes, link(), electronStore::links, linkTests_, and edm::OwnVector< T, P >::size().

Referenced by findBlocks().

218  {
219  constexpr unsigned rowsize = reco::PFBlockElement::kNBETypes;
220 
221  const edm::OwnVector<reco::PFBlockElement>& els = block.elements();
222 
223  block.bookLinkData();
224  unsigned elsize = els.size();
225  //First Loop: update all link data
226  for (unsigned i1 = 0; i1 < elsize; ++i1) {
227  for (unsigned i2 = 0; i2 < i1; ++i2) {
228  // no reflexive link
229  //if( i1==i2 ) continue;
230 
231  double dist = -1;
232 
233  bool linked = false;
234 
235  // are these elements already linked ?
236  // this can be optimized
237  const auto link_itr = links.find(std::make_pair(i2, i1));
238  if (link_itr != links.end()) {
239  dist = link_itr->second;
240  linked = true;
241  }
242 
243  if (!linked) {
244  const PFBlockElement::Type type1 = els[i1].type();
245  const PFBlockElement::Type type2 = els[i2].type();
246  const auto minmax = std::minmax(type1, type2);
247  const unsigned index = rowsize * minmax.second + minmax.first;
248  bool bTestLink =
249  (nullptr == linkTests_[index] ? false : linkTests_[index]->linkPrefilter(&(els[i1]), &(els[i2])));
250  if (bTestLink)
251  link(&els[i1], &els[i2], dist);
252  }
253 
254  //loading link data according to link test used: RECHIT
255  //block.setLink( i1, i2, chi2, block.linkData() );
256  block.setLink(i1, i2, dist, block.linkData());
257  }
258  }
259 }
size_type size() const
Definition: OwnVector.h:300
void link(const reco::PFBlockElement *el1, const reco::PFBlockElement *el2, double &dist) const
check whether 2 elements are linked. Returns distance
Definition: PFBlockAlgo.cc:261
std::vector< std::unique_ptr< BlockElementLinkerBase > > linkTests_
Definition: PFBlockAlgo.h:83

◆ setDebug()

void PFBlockAlgo::setDebug ( bool  debug)
inline

sets debug printout flag

Definition at line 59 of file PFBlockAlgo.h.

References debug, and debug_.

59 { debug_ = debug; }
#define debug
Definition: HDRShower.cc:19
bool debug_
if true, debug printouts activated
Definition: PFBlockAlgo.h:75

◆ setImporters()

void PFBlockAlgo::setImporters ( const std::vector< edm::ParameterSet > &  confs,
edm::ConsumesCollector cc 
)

Definition at line 124 of file PFBlockAlgo.cc.

References gpuPixelDoublets::cc, beamerCreator::create(), get, HLT_2024v10_cff::importerName, importers_, and AlCaHLTBitMon_QueryRunRegistry::string.

124  {
125  importers_.reserve(confs.size());
126  for (const auto& conf : confs) {
127  const std::string& importerName = conf.getParameter<std::string>("importerName");
129  }
130 }
def create(alignables, pedeDump, additionalData, outputFile, config)
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
std::vector< std::unique_ptr< BlockElementImporterBase > > importers_
Definition: PFBlockAlgo.h:80
#define get

◆ setLinkers()

void PFBlockAlgo::setLinkers ( const std::vector< edm::ParameterSet > &  confs)

Definition at line 84 of file PFBlockAlgo.cc.

References ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), beamerCreator::create(), elementTypes_, Exception, get, mps_fire::i, dqmiolumiharvest::j, kdtrees_, reco::PFBlockElement::kNBETypes, HLT_2024v10_cff::linkerName, linkTests_, linkTestSquare_, SiStripPI::max, SiStripPI::min, hcallasereventfilter2012_cfi::prefix, submitPVValidationJobs::split(), AlCaHLTBitMon_QueryRunRegistry::string, and HLT_2024v10_cff::useKDTree.

84  {
85  constexpr unsigned rowsize = reco::PFBlockElement::kNBETypes;
86  for (unsigned i = 0; i < rowsize; ++i) {
87  for (unsigned j = 0; j < rowsize; ++j) {
88  linkTestSquare_[i][j] = 0;
89  }
90  }
91  linkTests_.resize(rowsize * rowsize);
92  const std::string prefix("PFBlockElement::");
93  const std::string pfx_kdtree("KDTree");
94  for (const auto& conf : confs) {
95  const std::string& linkerName = conf.getParameter<std::string>("linkerName");
96  const std::string& linkTypeStr = conf.getParameter<std::string>("linkType");
97  size_t split = linkTypeStr.find(':');
98  if (split == std::string::npos) {
99  throw cms::Exception("MalformedLinkType") << "\"" << linkTypeStr << "\" is not a valid link type definition."
100  << " This string should have the form \"linkFrom:linkTo\"";
101  }
102  std::string link1(prefix + linkTypeStr.substr(0, split));
103  std::string link2(prefix + linkTypeStr.substr(split + 1, std::string::npos));
104  if (!(elementTypes_.count(link1) && elementTypes_.count(link2))) {
105  throw cms::Exception("InvalidBlockElementType")
106  << "One of \"" << link1 << "\" or \"" << link2 << "\" are invalid block element types!";
107  }
108  const PFBlockElement::Type type1 = elementTypes_.at(link1);
109  const PFBlockElement::Type type2 = elementTypes_.at(link2);
110  const unsigned index = rowsize * std::max(type1, type2) + std::min(type1, type2);
112  linkTestSquare_[type1][type2] = index;
113  linkTestSquare_[type2][type1] = index;
114  // setup KDtree if requested
115  const bool useKDTree = conf.getParameter<bool>("useKDTree");
116  if (useKDTree) {
117  kdtrees_.emplace_back(KDTreeLinkerFactory::get()->create(pfx_kdtree + linkerName, conf));
118  kdtrees_.back()->setTargetType(std::min(type1, type2));
119  kdtrees_.back()->setFieldType(std::max(type1, type2));
120  } // useKDTree
121  } // loop over confs
122 }
def create(alignables, pedeDump, additionalData, outputFile, config)
std::vector< std::unique_ptr< KDTreeLinkerBase > > kdtrees_
Definition: PFBlockAlgo.h:86
std::vector< std::unique_ptr< BlockElementLinkerBase > > linkTests_
Definition: PFBlockAlgo.h:83
const std::unordered_map< std::string, reco::PFBlockElement::Type > elementTypes_
Definition: PFBlockAlgo.h:82
unsigned int linkTestSquare_[reco::PFBlockElement::kNBETypes][reco::PFBlockElement::kNBETypes]
Definition: PFBlockAlgo.h:84
#define get

◆ updateEventSetup()

void PFBlockAlgo::updateEventSetup ( const edm::EventSetup es)

Definition at line 275 of file PFBlockAlgo.cc.

References importers_.

Referenced by PFBlockProducer::beginLuminosityBlock().

275  {
276  for (auto& importer : importers_) {
277  importer->updateEventSetup(es);
278  }
279 }
std::vector< std::unique_ptr< BlockElementImporterBase > > importers_
Definition: PFBlockAlgo.h:80

Friends And Related Function Documentation

◆ operator<<

std::ostream& operator<< ( std::ostream &  out,
const PFBlockAlgo a 
)
friend

Definition at line 325 of file PFBlockAlgo.cc.

325  {
326  if (!out)
327  return out;
328 
329  out << "====== Particle Flow Block Algorithm ======= ";
330  out << endl;
331  out << "number of unassociated elements : " << a.elements_.size() << endl;
332  out << endl;
333 
334  for (auto const& element : a.elements_) {
335  out << "\t" << *element << endl;
336  }
337 
338  return out;
339 }
double a
Definition: hdecay.h:121

Member Data Documentation

◆ debug_

bool PFBlockAlgo::debug_
private

if true, debug printouts activated

Definition at line 75 of file PFBlockAlgo.h.

Referenced by link(), setDebug(), and ~PFBlockAlgo().

◆ elements_

ElementList PFBlockAlgo::elements_
private

Definition at line 71 of file PFBlockAlgo.h.

Referenced by buildElements(), findBlocks(), and ~PFBlockAlgo().

◆ elementTypes_

const std::unordered_map<std::string, reco::PFBlockElement::Type> PFBlockAlgo::elementTypes_
private

Definition at line 82 of file PFBlockAlgo.h.

Referenced by setLinkers().

◆ importers_

std::vector<std::unique_ptr<BlockElementImporterBase> > PFBlockAlgo::importers_
private

Definition at line 80 of file PFBlockAlgo.h.

Referenced by buildElements(), setImporters(), and updateEventSetup().

◆ kdtrees_

std::vector<std::unique_ptr<KDTreeLinkerBase> > PFBlockAlgo::kdtrees_
private

Definition at line 86 of file PFBlockAlgo.h.

Referenced by buildElements(), findBlocks(), and setLinkers().

◆ linkTests_

std::vector<std::unique_ptr<BlockElementLinkerBase> > PFBlockAlgo::linkTests_
private

Definition at line 83 of file PFBlockAlgo.h.

Referenced by findBlocks(), link(), packLinks(), and setLinkers().

◆ linkTestSquare_

unsigned int PFBlockAlgo::linkTestSquare_[reco::PFBlockElement::kNBETypes][reco::PFBlockElement::kNBETypes]
private

Definition at line 84 of file PFBlockAlgo.h.

Referenced by findBlocks(), and setLinkers().

◆ ranges_

ElementRanges PFBlockAlgo::ranges_
private

Definition at line 72 of file PFBlockAlgo.h.

Referenced by buildElements(), and findBlocks().

◆ useHO_

bool PFBlockAlgo::useHO_
private

Definition at line 78 of file PFBlockAlgo.h.