CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFBlockAlgo.cc
Go to the documentation of this file.
7 
9 
10 #include <stdexcept>
11 #include <algorithm>
12 #include "TMath.h"
13 
14 using namespace std;
15 using namespace reco;
16 
17 #define INIT_ENTRY(name) {#name,name}
18 
19 namespace {
20  class QuickUnion{
21  std::vector<unsigned> _id;
22  std::vector<unsigned> _size;
23  int _count;
24 
25  public:
26  QuickUnion(const unsigned NBranches) {
27  _count = NBranches;
28  _id.resize(NBranches);
29  _size.resize(NBranches);
30  for( unsigned i = 0; i < NBranches; ++i ) {
31  _id[i] = i;
32  _size[i] = 1;
33  }
34  }
35 
36  int count() const { return _count; }
37 
38  unsigned find(unsigned p) {
39  while( p != _id[p] ) {
40  _id[p] = _id[_id[p]];
41  p = _id[p];
42  }
43  return p;
44  }
45 
46  bool connected(unsigned p, unsigned q) { return find(p) == find(q); }
47 
48  void unite(unsigned p, unsigned q) {
49  unsigned rootP = find(p);
50  unsigned rootQ = find(q);
51  _id[p] = q;
52 
53  if(_size[rootP] < _size[rootQ] ) {
54  _id[rootP] = rootQ; _size[rootQ] += _size[rootP];
55  } else {
56  _id[rootQ] = rootP; _size[rootP] += _size[rootQ];
57  }
58  --_count;
59  }
60  };
61 }
62 
63 
64 //for debug only
65 //#define PFLOW_DEBUG
66 
68  blocks_( new reco::PFBlockCollection ),
69  debug_(false),
70  _elementTypes( {
71  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  } ) {}
83 
84 void PFBlockAlgo::setLinkers(const std::vector<edm::ParameterSet>& confs) {
85  constexpr unsigned rowsize = reco::PFBlockElement::kNBETypes;
86  for( unsigned i = 0; i < rowsize; ++i ) {
87  for( unsigned j = 0; j < rowsize; ++j ) {
88 
89  _linkTestSquare[i][j] = 0;
90  }
91  }
92  _linkTests.resize(rowsize*rowsize);
93  const std::string prefix("PFBlockElement::");
94  const std::string pfx_kdtree("KDTree");
95  for( const auto& conf : confs ) {
96  const std::string& linkerName =
97  conf.getParameter<std::string>("linkerName");
98  const std::string& linkTypeStr =
99  conf.getParameter<std::string>("linkType");
100  size_t split = linkTypeStr.find(':');
101  if( split == std::string::npos ) {
102  throw cms::Exception("MalformedLinkType")
103  << "\"" << linkTypeStr << "\" is not a valid link type definition."
104  << " This string should have the form \"linkFrom:linkTo\"";
105  }
106  std::string link1(prefix+linkTypeStr.substr(0,split));
107  std::string link2(prefix+linkTypeStr.substr(split+1,std::string::npos));
108  if( !(_elementTypes.count(link1) && _elementTypes.count(link2) ) ) {
109  throw cms::Exception("InvalidBlockElementType")
110  << "One of \"" << link1 << "\" or \"" << link2
111  << "\" are invalid block element types!";
112  }
113  const PFBlockElement::Type type1 = _elementTypes.at(link1);
114  const PFBlockElement::Type type2 = _elementTypes.at(link2);
115  const unsigned index = rowsize*std::max(type1,type2)+std::min(type1,type2);
116  BlockElementLinkerBase * linker =
117  BlockElementLinkerFactory::get()->create(linkerName,conf);
118  _linkTests[index].reset(linker);
119  _linkTestSquare[type1][type2] = index;
120  _linkTestSquare[type2][type1] = index;
121  // setup KDtree if requested
122  const bool useKDTree = conf.getParameter<bool>("useKDTree");
123  if( useKDTree ) {
124  _kdtrees.emplace_back( KDTreeLinkerFactory::get()->create(pfx_kdtree+
125  linkerName) );
126  _kdtrees.back()->setTargetType(std::min(type1,type2));
127  _kdtrees.back()->setFieldType(std::max(type1,type2));
128  }
129  }
130 }
131 
132 void PFBlockAlgo::setImporters(const std::vector<edm::ParameterSet>& confs,
133  edm::ConsumesCollector& sumes) {
134  _importers.reserve(confs.size());
135  for( const auto& conf : confs ) {
136  const std::string& importerName =
137  conf.getParameter<std::string>("importerName");
138  BlockElementImporterBase * importer =
139  BlockElementImporterFactory::get()->create(importerName,conf,sumes);
140  _importers.emplace_back(importer);
141  }
142 }
143 
145 
146 #ifdef PFLOW_DEBUG
147  if(debug_)
148  cout<<"~PFBlockAlgo - number of remaining elements: "
149  <<elements_.size()<<endl;
150 #endif
151 }
152 
154  // Glowinski & Gouzevitch
155  for( const auto& kdtree : _kdtrees ) {
156  kdtree->process();
157  }
158  // !Glowinski & Gouzevitch
159  // the blocks have not been passed to the event, and need to be cleared
160  if( blocks_.get() ) blocks_->clear();
161  else blocks_.reset( new reco::PFBlockCollection );
162  blocks_->reserve(elements_.size());
163 
164  QuickUnion qu(bare_elements_.size());
165  const auto elem_size = bare_elements_.size();
166  for( unsigned i = 0; i < elem_size; ++i ) {
167  for( unsigned j = 0; j < elem_size; ++j ) {
168  if( qu.connected(i,j) || j == i ) continue;
170  j = ranges_[bare_elements_[j]->type()].second;
171  continue;
172  }
173  auto p1(bare_elements_[i]), p2(bare_elements_[j]);
174  const PFBlockElement::Type type1 = p1->type();
175  const PFBlockElement::Type type2 = p2->type();
176  const unsigned index = _linkTestSquare[type1][type2];
177  if( _linkTests[index]->linkPrefilter(p1,p2) ) {
178  const double dist = _linkTests[index]->testLink(p1,p2);
179  // compute linking info if it is possible
180  if( dist > -0.5 ) {
181  qu.unite(i,j);
182  }
183  }
184  }
185  }
186 
187  std::unordered_multimap<unsigned,unsigned> blocksmap(elements_.size());
188  std::vector<unsigned> keys;
189  keys.reserve(elements_.size());
190  for( unsigned i = 0; i < elements_.size(); ++i ) {
191  unsigned key = i;
192  while( key != qu.find(key) ) key = qu.find(key); // make sure we always find the root node...
193  auto pos = std::lower_bound(keys.begin(),keys.end(),key);
194  if( pos == keys.end() || *pos != key ) {
195  keys.insert(pos,key);
196  }
197  blocksmap.emplace(key,i);
198  }
199 
201  PFBlock::LinkTest linktest = PFBlock::LINKTEST_RECHIT;
202  for( auto key : keys ) {
203  blocks_->push_back( reco::PFBlock() );
204  auto range = blocksmap.equal_range(key);
205  auto& the_block = blocks_->back();
206  ElementList::value_type::pointer p1(bare_elements_[range.first->second]);
207  the_block.addElement(p1);
208  const unsigned block_size = blocksmap.count(key) + 1;
209  //reserve up to 1M or 8MB; pay rehash cost for more
210  std::unordered_map<std::pair<unsigned int,unsigned int>, PFBlockLink > links(min(1000000u,block_size*block_size));
211  auto itr = range.first;
212  ++itr;
213  for( ; itr != range.second; ++itr ) {
214  ElementList::value_type::pointer p2(bare_elements_[itr->second]);
215  const PFBlockElement::Type type1 = p1->type();
216  const PFBlockElement::Type type2 = p2->type();
217  the_block.addElement(p2);
218  linktest = PFBlock::LINKTEST_RECHIT; //rechit by default
219  linktype = static_cast<PFBlockLink::Type>(1<<(type1-1)|1<<(type2-1));
220  const unsigned index = _linkTestSquare[type1][type2];
221  if( nullptr != _linkTests[index] ) {
222  const double dist = _linkTests[index]->testLink(p1,p2);
223  links.emplace( std::make_pair(p1->index(), p2->index()) ,
224  PFBlockLink( linktype, linktest, dist,
225  p1->index(), p2->index() ) );
226  }
227  }
228  packLinks( the_block, links );
229  }
230 
231  bare_elements_.clear();
232  elements_.clear();
233 }
234 
235 void
237  const std::unordered_map<std::pair<unsigned int,unsigned int>,PFBlockLink>& links ) const {
238  constexpr unsigned rowsize = reco::PFBlockElement::kNBETypes;
239 
241 
242  block.bookLinkData();
243  unsigned elsize = els.size();
244  //First Loop: update all link data
245  for( unsigned i1=0; i1<elsize; ++i1 ) {
246  for( unsigned i2=0; i2<i1; ++i2 ) {
247 
248  // no reflexive link
249  //if( i1==i2 ) continue;
250 
251  double dist = -1;
252 
253  bool linked = false;
254  PFBlock::LinkTest linktest
255  = PFBlock::LINKTEST_RECHIT;
256 
257  // are these elements already linked ?
258  // this can be optimized
259  const auto link_itr = links.find(std::make_pair(i2,i1));
260  if( link_itr != links.end() ) {
261  dist = link_itr->second.dist();
262  linktest = link_itr->second.test();
263  linked = true;
264  }
265 
266  if(!linked) {
267  const PFBlockElement::Type type1 = els[i1].type();
268  const PFBlockElement::Type type2 = els[i2].type();
269  const auto minmax = std::minmax(type1,type2);
270  const unsigned index = rowsize*minmax.second + minmax.first;
272  bool bTestLink = ( nullptr == _linkTests[index] ? false : _linkTests[index]->linkPrefilter(&(els[i1]),&(els[i2])) );
273  if (bTestLink) link( & els[i1], & els[i2], linktype, linktest, dist);
274  }
275 
276  //loading link data according to link test used: RECHIT
277  //block.setLink( i1, i2, chi2, block.linkData() );
278 #ifdef PFLOW_DEBUG
279  if( debug_ )
280  cout << "Setting link between elements " << i1 << " and " << i2
281  << " of dist =" << dist << " computed from link test "
282  << linktest << endl;
283 #endif
284  block.setLink( i1, i2, dist, block.linkData(), linktest );
285  }
286  }
287 
288 }
289 
290 // see plugins/linkers for the functions that calculate distances
291 // for each available link type
292 inline bool
294  const reco::PFBlockElement* next) const {
295  constexpr unsigned rowsize = reco::PFBlockElement::kNBETypes;
296  const PFBlockElement::Type type1 = (last)->type();
297  const PFBlockElement::Type type2 = (next)->type();
298  const unsigned index = rowsize*std::max(type1,type2) + std::min(type1,type2);
299  bool result = _linkTests[index]->linkPrefilter(last,next);
300  return result;
301 }
302 
303 inline void
305  const reco::PFBlockElement* el2,
306  PFBlockLink::Type& linktype,
307  reco::PFBlock::LinkTest& linktest,
308  double& dist) const {
309  constexpr unsigned rowsize = reco::PFBlockElement::kNBETypes;
310  dist=-1.0;
311  linktest = PFBlock::LINKTEST_RECHIT; //rechit by default
312  const PFBlockElement::Type type1 = el1->type();
313  const PFBlockElement::Type type2 = el2->type();
314  linktype = static_cast<PFBlockLink::Type>(1<<(type1-1)|1<<(type2-1));
315  const unsigned index = rowsize*std::max(type1,type2) + std::min(type1,type2);
316  if(debug_ ) {
317  std::cout << " PFBlockAlgo links type1 " << type1
318  << " type2 " << type2 << std::endl;
319  }
320 
321  // index is always checked in the preFilter above, no need to check here
322  dist = _linkTests[index]->testLink(el1,el2);
323 }
324 
326  for( auto& importer : _importers ) {
327  importer->updateEventSetup(es);
328  }
329 }
330 
331 // see plugins/importers and plugins/kdtrees
332 // for the definitions of available block element importers
333 // and kdtree preprocessors
334 void PFBlockAlgo::buildElements(const edm::Event& evt) {
335  // import block elements as defined in python configuration
336  ranges_.fill(std::make_pair(0,0));
337  elements_.clear();
338  for( const auto& importer : _importers ) {
339  importer->importToBlock(evt,elements_);
340  }
341 
342  std::sort(elements_.begin(),elements_.end(),
343  [](const auto& a, const auto& b) { return a->type() < b->type(); } );
344 
345  bare_elements_.resize(elements_.size());
346  for( unsigned i = 0; i < elements_.size(); ++i ) {
347  bare_elements_[i] = elements_[i].get();
348  }
349 
350  // list is now partitioned, so mark the boundaries so we can efficiently skip chunks
351  unsigned current_type = ( elements_.size() ? elements_[0]->type() : 0 );
352  unsigned last_type = ( elements_.size() ? elements_.back()->type() : 0 );
353  ranges_[current_type].first = 0;
354  ranges_[last_type].second = elements_.size()-1;
355  for( size_t i = 0; i < elements_.size(); ++i ) {
356  const auto the_type = elements_[i]->type();
357  if( the_type != current_type ) {
358  ranges_[the_type].first = i;
359  ranges_[current_type].second = i-1;
360  current_type = the_type;
361  }
362  }
363  // -------------- Loop over block elements ---------------------
364 
365  // Here we provide to all KDTree linkers the collections to link.
366  // Glowinski & Gouzevitch
367 
368  for (ElementList::iterator it = elements_.begin();
369  it != elements_.end(); ++it) {
370  for( const auto& kdtree : _kdtrees ) {
371  if( (*it)->type() == kdtree->targetType() ) {
372  kdtree->insertTargetElt(it->get());
373  }
374  if( (*it)->type() == kdtree->fieldType() ) {
375  kdtree->insertFieldClusterElt(it->get());
376  }
377  }
378  }
379  //std::cout << "(new) imported: " << elements_.size() << " elements!" << std::endl;
380 }
381 
382 std::ostream& operator<<(std::ostream& out, const PFBlockAlgo& a) {
383  if(! out) return out;
384 
385  out<<"====== Particle Flow Block Algorithm ======= ";
386  out<<endl;
387  out<<"number of unassociated elements : "<<a.elements_.size()<<endl;
388  out<<endl;
389 
390  for(PFBlockAlgo::IEC ie = a.elements_.begin();
391  ie != a.elements_.end(); ++ie) {
392  out<<"\t"<<**ie <<endl;
393  }
394 
395 
396  // const PFBlockCollection& blocks = a.blocks();
397 
398  const std::auto_ptr< reco::PFBlockCollection >& blocks
399  = a.blocks();
400 
401  if(!blocks.get() ) {
402  out<<"blocks already transfered"<<endl;
403  }
404  else {
405  out<<"number of blocks : "<<blocks->size()<<endl;
406  out<<endl;
407 
408  for(PFBlockAlgo::IBC ib=blocks->begin();
409  ib != blocks->end(); ++ib) {
410  out<<(*ib)<<endl;
411  }
412  }
413 
414  return out;
415 }
416 
417 // a little history, ideas we may want to keep around for later
418  /*
419  // Links between the two preshower layers are not used for now - disable
420  case PFBlockLink::PS1andPS2:
421  {
422  PFClusterRef ps1ref = lowEl->clusterRef();
423  PFClusterRef ps2ref = highEl->clusterRef();
424  assert( !ps1ref.isNull() );
425  assert( !ps2ref.isNull() );
426  // PJ - 14-May-09 : A link by rechit is needed here !
427  // dist = testPS1AndPS2( *ps1ref, *ps2ref );
428  dist = -1.;
429  linktest = PFBlock::LINKTEST_RECHIT;
430  break;
431  }
432  case PFBlockLink::TRACKandPS1:
433  case PFBlockLink::TRACKandPS2:
434  {
435  //cout<<"TRACKandPS"<<endl;
436  PFRecTrackRef trackref = lowEl->trackRefPF();
437  PFClusterRef clusterref = highEl->clusterRef();
438  assert( !trackref.isNull() );
439  assert( !clusterref.isNull() );
440  // PJ - 14-May-09 : A link by rechit is needed here !
441  dist = testTrackAndPS( *trackref, *clusterref );
442  linktest = PFBlock::LINKTEST_RECHIT;
443  break;
444  }
445  // GSF Track/Brem Track and preshower cluster links are not used for now - disable
446  case PFBlockLink::PS1andGSF:
447  case PFBlockLink::PS2andGSF:
448  {
449  PFClusterRef psref = lowEl->clusterRef();
450  assert( !psref.isNull() );
451  const reco::PFBlockElementGsfTrack * GsfEl = dynamic_cast<const reco::PFBlockElementGsfTrack*>(highEl);
452  const PFRecTrack * myTrack = &(GsfEl->GsftrackPF());
453  // PJ - 14-May-09 : A link by rechit is needed here !
454  dist = testTrackAndPS( *myTrack, *psref );
455  linktest = PFBlock::LINKTEST_RECHIT;
456  break;
457  }
458  case PFBlockLink::PS1andBREM:
459  case PFBlockLink::PS2andBREM:
460  {
461  PFClusterRef psref = lowEl->clusterRef();
462  assert( !psref.isNull() );
463  const reco::PFBlockElementBrem * BremEl = dynamic_cast<const reco::PFBlockElementBrem*>(highEl);
464  const PFRecTrack * myTrack = &(BremEl->trackPF());
465  // PJ - 14-May-09 : A link by rechit is needed here !
466  dist = testTrackAndPS( *myTrack, *psref );
467  linktest = PFBlock::LINKTEST_RECHIT;
468  break;
469  }
470  */
const std::auto_ptr< reco::PFBlockCollection > & blocks() const
Definition: PFBlockAlgo.h:129
type
Definition: HCALResponse.h:21
Abstract base class for a PFBlock element (track, cluster...)
int i
Definition: DBlmapReader.cc:9
std::vector< LinkTestPtr > _linkTests
Definition: PFBlockAlgo.h:173
const std::unordered_map< std::string, reco::PFBlockElement::Type > _elementTypes
Definition: PFBlockAlgo.h:172
int ib
Definition: cuy.py:660
ElementList::const_iterator IEC
Definition: PFBlockAlgo.h:101
Type type() const
size_type size() const
Definition: OwnVector.h:264
const edm::OwnVector< reco::PFBlockElement > & elements() const
Definition: PFBlock.h:107
std::auto_ptr< reco::PFBlockCollection > blocks_
Definition: PFBlockAlgo.h:156
const LinkData & linkData() const
Definition: PFBlock.h:112
std::ostream & operator<<(std::ostream &out, const ALILine &li)
Definition: ALILine.cc:188
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
ElementList elements_
Definition: PFBlockAlgo.h:159
#define constexpr
void setLink(unsigned i1, unsigned i2, double dist, LinkData &linkData, LinkTest test=LINKTEST_RECHIT) const
Definition: PFBlock.cc:26
tuple result
Definition: mps_fire.py:84
Particle Flow Algorithm.
Definition: PFBlockAlgo.h:91
#define INIT_ENTRY(name)
Definition: PFBlockAlgo.cc:17
void link(const reco::PFBlockElement *el1, const reco::PFBlockElement *el2, PFBlockLink::Type &linktype, reco::PFBlock::LinkTest &linktest, double &dist) const
check whether 2 elements are linked. Returns distance and linktype
std::vector< KDTreePtr > _kdtrees
Definition: PFBlockAlgo.h:176
int j
Definition: DBlmapReader.cc:9
void updateEventSetup(const edm::EventSetup &)
void setImporters(const std::vector< edm::ParameterSet > &, edm::ConsumesCollector &)
T min(T a, T b)
Definition: MathUtil.h:58
std::vector< PFBlock > PFBlockCollection
collection of PFBlock objects
Definition: PFBlockFwd.h:11
void bookLinkData()
Definition: PFBlock.cc:20
double p2[4]
Definition: TauolaWrapper.h:90
string key
FastSim: produces sample of signal events, overlayed with premixed minbias events.
ElementRanges ranges_
Definition: PFBlockAlgo.h:161
void packLinks(reco::PFBlock &block, const std::unordered_map< std::pair< unsigned int, unsigned int >, PFBlockLink > &links) const
reco::PFBlockCollection::const_iterator IBC
Definition: PFBlockAlgo.h:102
list blocks
Definition: gather_cfg.py:90
double b
Definition: hdecay.h:120
bool linkPrefilter(const reco::PFBlockElement *last, const reco::PFBlockElement *next) const
Avoid to check links when not useful.
void findBlocks()
build blocks
double p1[4]
Definition: TauolaWrapper.h:89
double a
Definition: hdecay.h:121
std::vector< ElementList::value_type::pointer > bare_elements_
Definition: PFBlockAlgo.h:160
tuple cout
Definition: gather_cfg.py:145
volatile std::atomic< bool > shutdown_flag false
void buildElements(const edm::Event &)
unsigned int _linkTestSquare[reco::PFBlockElement::kNBETypes][reco::PFBlockElement::kNBETypes]
Definition: PFBlockAlgo.h:174
bool debug_
if true, debug printouts activated
Definition: PFBlockAlgo.h:164
void setLinkers(const std::vector< edm::ParameterSet > &)
double split
Definition: MVATrainer.cc:139
T get(const Candidate &c)
Definition: component.h:55
std::vector< ImporterPtr > _importers
Definition: PFBlockAlgo.h:169
Block of elements.
Definition: PFBlock.h:30