CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
clusteringTreeConverters.h
Go to the documentation of this file.
1 //
2 // Code converting data between FFTJet clustering tree objects and
3 // event-storable entities. The data can be stored in either single
4 // or double precision. If the data is stored in double precision,
5 // the clustering tree can be restored in its original form.
6 //
7 // Written by:
8 //
9 // Terence Libeiro
10 // Igor Volobouev
11 //
12 // June 2010
13 
14 #ifndef RecoJets_FFTJetAlgorithms_clusteringTreeConverters_h
15 #define RecoJets_FFTJetAlgorithms_clusteringTreeConverters_h
16 
17 #include <cfloat>
18 
19 #include "fftjet/Peak.hh"
20 #include "fftjet/AbsClusteringTree.hh"
21 #include "fftjet/SparseClusteringTree.hh"
22 
26 
27 namespace fftjetcms {
28  // The function below makes PattRecoTree storable in the event
29  // record out of a sparse clustering tree of fftjet::Peak objects
30  template<class Real>
32  const fftjet::SparseClusteringTree<fftjet::Peak,long>& in,
33  bool writeOutScaleInfo,
35 
36  // The function below restores sparse clustering tree of fftjet::Peak
37  // objects using the PattRecoTree. The "completeEventScale" parameter
38  // should be set to the appropriate scale in case the complete event
39  // was written out and its scale not stored in the event. If the complete
40  // event was not written out, the scale must be set to 0.
41  template<class Real>
44  const std::vector<double>* scaleSetIfNotAdaptive,
45  double completeEventScale,
46  fftjet::SparseClusteringTree<fftjet::Peak,long>* out);
47 
48  // The function below makes PattRecoTree storable in the event
49  // record out of a dense clustering tree of fftjet::Peak objects
50  template<class Real>
52  const fftjet::AbsClusteringTree<fftjet::Peak,long>& in,
53  bool writeOutScaleInfo,
55 
56  // The function below restores dense clustering tree of fftjet::Peak
57  // objects using the PattRecoTree. The "completeEventScale" parameter
58  // should be set to the appropriate scale in case the complete event
59  // was written out and its scale not stored in the event. If the complete
60  // event was not written out, the scale must be set to 0.
61  template<class Real>
64  const std::vector<double>* scaleSetIfNotAdaptive,
65  double completeEventScale,
66  fftjet::AbsClusteringTree<fftjet::Peak,long>* out);
67 }
68 
70 //
71 // Implementation follows
72 //
74 
75 namespace fftjetcms {
76  // The function below makes PattRecoTree storable in the event
77  // record out of a sparse clustering tree of fftjet::Peak objects
78  template<class Real>
80  const fftjet::SparseClusteringTree<fftjet::Peak,long>& sparseTree,
81  const bool writeOutScaleInfo,
83  {
84  typedef fftjet::SparseClusteringTree<fftjet::Peak,long> SparseTree;
85  typedef reco::PattRecoPeak<Real> StoredPeak;
86  typedef reco::PattRecoNode<StoredPeak> StoredNode;
87  typedef reco::PattRecoTree<Real,StoredPeak> StoredTree;
88 
89  assert(tree);
90 
91  tree->clear();
92  tree->setSparse(true);
93 
94  const unsigned nNodes = sparseTree.size();
95  double hessian[3] = {0., 0., 0.};
96 
97  // Do not write out the meaningless top node
98  tree->reserveNodes(nNodes - 1);
99  for (unsigned i=1; i<nNodes; ++i)
100  {
101  const SparseTree::Node& node(sparseTree.getNode(i));
102  const fftjet::Peak& peak(node.getCluster());
103  peak.hessian(hessian);
104  StoredNode sn(StoredPeak(peak.eta(),
105  peak.phi(),
106  peak.magnitude(),
107  hessian,
108  peak.driftSpeed(),
109  peak.magSpeed(),
110  peak.lifetime(),
111  peak.scale(),
112  peak.nearestNeighborDistance(),
113  peak.clusterRadius(),
114  peak.clusterSeparation()),
115  node.originalLevel(),
116  node.mask(),
117  node.parent());
118  tree->addNode(sn);
119  }
120 
121  // Do we want to write out the scales? We will use the following
122  // convention: if the tree is using an adaptive algorithm, the scales
123  // will be written out. If not, they are not going to change from
124  // event to event. In this case the scales would waste disk space
125  // for no particularly good reason, so they will not be written out.
126  if (writeOutScaleInfo)
127  {
128  // Do not write out the meaningless top-level scale
129  const unsigned nScales = sparseTree.nLevels();
130  tree->reserveScales(nScales - 1);
131  for (unsigned i=1; i<nScales; ++i)
132  tree->addScale(sparseTree.getScale(i));
133  }
134  }
135 
136  // The function below restores sparse clustering tree of fftjet::Peak
137  // objects using the PattRecoTree
138  template<class Real>
141  const std::vector<double>* scaleSetIfNotAdaptive,
142  const double completeEventScale,
143  fftjet::SparseClusteringTree<fftjet::Peak,long>* out)
144  {
145  typedef fftjet::SparseClusteringTree<fftjet::Peak,long> SparseTree;
146  typedef reco::PattRecoPeak<Real> StoredPeak;
147  typedef reco::PattRecoNode<StoredPeak> StoredNode;
148  typedef reco::PattRecoTree<Real,StoredPeak> StoredTree;
149 
150  if (!in.isSparse())
151  throw cms::Exception("FFTJetBadConfig")
152  << "can't restore sparse clustering tree"
153  << " from densely stored record" << std::endl;
154 
155  assert(out);
156  out->clear();
157 
158  const std::vector<StoredNode>& nodes(in.getNodes());
159  const unsigned n = nodes.size();
160  out->reserveNodes(n + 1U);
161 
162  double hessian[3] = {0., 0., 0.};
163 
164  for (unsigned i=0; i<n; ++i)
165  {
166  const StoredNode& snode(nodes[i]);
167  const StoredPeak& p(snode.getCluster());
168  p.hessian(hessian);
169  const SparseTree::Node node(
170  fftjet::Peak(p.eta(), p.phi(), p.magnitude(),
171  hessian, p.driftSpeed(),
172  p.magSpeed(), p.lifetime(),
173  p.scale(), p.nearestNeighborDistance(),
174  1.0, 0.0, 0.0,
175  p.clusterRadius(), p.clusterSeparation()),
176  snode.originalLevel(),
177  snode.mask());
178  out->addNode(node, snode.parent());
179  }
180 
181  const std::vector<Real>& storedScales(in.getScales());
182  if (!storedScales.empty())
183  {
184  const unsigned nsc = storedScales.size();
185  out->reserveScales(nsc + 1U);
186  out->addScale(DBL_MAX);
187  const Real* scales = &storedScales[0];
188  for (unsigned i=0; i<nsc; ++i)
189  out->addScale(scales[i]);
190  }
191  else if (scaleSetIfNotAdaptive && !scaleSetIfNotAdaptive->empty())
192  {
193  const unsigned nsc = scaleSetIfNotAdaptive->size();
194  // There may be the "complete event" scale added at the end.
195  // Reserve a sufficient number of scales to take this into
196  // account.
197  if (completeEventScale)
198  out->reserveScales(nsc + 2U);
199  else
200  out->reserveScales(nsc + 1U);
201  out->addScale(DBL_MAX);
202  const double* scales = &(*scaleSetIfNotAdaptive)[0];
203  for (unsigned i=0; i<nsc; ++i)
204  out->addScale(scales[i]);
205  if (completeEventScale)
206  out->addScale(completeEventScale);
207  }
208  else
209  {
210  throw cms::Exception("FFTJetBadConfig")
211  << "can't restore sparse clustering tree scales"
212  << std::endl;
213  }
214  }
215 
216  // The function below makes PattRecoTree storable in the event
217  // record out of a dense clustering tree of fftjet::Peak objects
218  template<class Real>
220  const fftjet::AbsClusteringTree<fftjet::Peak,long>& in,
221  bool writeOutScaleInfo,
223  {
224  typedef fftjet::AbsClusteringTree<fftjet::Peak,long> DenseTree;
225  typedef reco::PattRecoPeak<Real> StoredPeak;
226  typedef reco::PattRecoNode<StoredPeak> StoredNode;
227  typedef reco::PattRecoTree<Real,StoredPeak> StoredTree;
228 
229 
230  assert(out);
231  out->clear();
232  out->setSparse(false);
233 
234 
235  const unsigned nLevels = in.nLevels();
236  double hessian[3] = {0., 0., 0.};
237 
238 
239  // Do not write out the meaningless top node
240  out->reserveNodes(in.nClusters() - 1);
241 
242 
243  for (unsigned i=1; i<nLevels; ++i)
244  {
245 
246  const unsigned int nclus = in.nClusters(i);
247  DenseTree::NodeId id(i,0);
248  for (;id.second<nclus; ++id.second)
249  {
250 
251 
252  const fftjet::Peak& peak(in.getCluster(id));
253  peak.hessian(hessian);
254  StoredNode sn(StoredPeak(peak.eta(),
255  peak.phi(),
256  peak.magnitude(),
257  hessian,
258  peak.driftSpeed(),
259  peak.magSpeed(),
260  peak.lifetime(),
261  peak.scale(),
262  peak.nearestNeighborDistance(),
263  peak.clusterRadius(),
264  peak.clusterSeparation()),
265  i,
266  0,
267  0);
268 
269  out->addNode(sn);
270 
271  }
272 
273  }
274 
275 
276  // Do we want to write out the scales? We will use the following
277  // convention: if the tree is using an adaptive algorithm, the scales
278  // will be written out. If not, they are not going to change from
279  // event to event. In this case the scales would waste disk space
280  // for no particularly good reason, so they will not be written out.
281  if (writeOutScaleInfo)
282  {
283  // Do not write out the meaningless top-level scale
284  const unsigned nScales = in.nLevels();
285 
286  out->reserveScales(nScales - 1);
287 
288  for (unsigned i=1; i<nScales; ++i)
289  out->addScale(in.getScale(i));
290 
291  }
292  }
293 
294  // The function below restores dense clustering tree of fftjet::Peak
295  // objects using the PattRecoTree
296  template<class Real>
299  const std::vector<double>* scaleSetIfNotAdaptive,
300  const double completeEventScale,
301  fftjet::AbsClusteringTree<fftjet::Peak,long>* out)
302  {
303  typedef fftjet::AbsClusteringTree<fftjet::Peak,long> DenseTree;
304  typedef reco::PattRecoPeak<Real> StoredPeak;
305  typedef reco::PattRecoNode<StoredPeak> StoredNode;
306  typedef reco::PattRecoTree<Real,StoredPeak> StoredTree;
307 
308  if (in.isSparse())
309  throw cms::Exception("FFTJetBadConfig")
310  << "can't restore dense clustering tree"
311  << " from sparsely stored record" << std::endl;
312 
313  assert(out);
314  out->clear();
315 
316  const std::vector<StoredNode>& nodes(in.getNodes());
317  const unsigned n = nodes.size();
318  double hessian[3] = {0., 0., 0.};
319 
320  const std::vector<Real>& scales (in.getScales());
321  unsigned int scnum = 0;
322  std::vector<fftjet::Peak> clusters;
323  const unsigned scsize1 = scales.size();
324  unsigned scsize2 = scaleSetIfNotAdaptive ? scaleSetIfNotAdaptive->size() : 0;
325  if (scsize2 && completeEventScale) ++scsize2;
326  const unsigned scsize = (scsize1==0?scsize2:scsize1);
327 
328  if (scsize == 0)
329  throw cms::Exception("FFTJetBadConfig")
330  << " No scales passed to the function densePeakTreeFromStorable()"
331  << std::endl;
332 
333  // to check whether the largest level equals the size of scale vector
334  const double* sc_not_ad = scsize2 ? &(*scaleSetIfNotAdaptive)[0] : 0;
335 
336  unsigned templevel = n ? nodes[0].originalLevel() : 1;
337  for (unsigned i=0; i<n; ++i)
338  {
339  const StoredNode& snode(nodes[i]);
340  const StoredPeak& p(snode.getCluster());
341  p.hessian(hessian);
342 
343  const unsigned levelNumber = snode.originalLevel();
344 
345  if (templevel != levelNumber)
346  {
347  if (scnum >= scsize)
348  throw cms::Exception("FFTJetBadConfig")
349  << "bad scales, please check the scales"
350  << std::endl;
351  const double scale = ( (scsize1==0) ? sc_not_ad[scnum] : scales[scnum] );
352  out->insert(scale, clusters, 0L);
353  clusters.clear();
354  templevel = levelNumber;
355  ++scnum;
356  }
357 
358  fftjet::Peak apeak(p.eta(), p.phi(), p.magnitude(),
359  hessian, p.driftSpeed(),
360  p.magSpeed(), p.lifetime(),
361  p.scale(), p.nearestNeighborDistance(),
362  1.0, 0.0, 0.0,
363  p.clusterRadius(), p.clusterSeparation());
364  clusters.push_back(apeak);
365 
366  if (i==(n-1) && levelNumber!=scsize)
367  throw cms::Exception("FFTJetBadConfig")
368  << "bad scales, please check the scales"
369  << std::endl;
370  }
371 
372  const double scale = scsize1 ? scales[scnum] : completeEventScale ?
373  completeEventScale : sc_not_ad[scnum];
374  out->insert(scale, clusters, 0L);
375  }
376 }
377 
378 #endif // RecoJets_FFTJetAlgorithms_clusteringTreeConverters_h
int i
Definition: DBlmapReader.cc:9
Tree nodes for storing FFTJet preclusters.
Definition: PattRecoNode.h:16
Class for storing FFTJet sparse clustering trees.
Definition: PattRecoTree.h:18
Preclusters from FFTJet pattern recognition stage.
Definition: PattRecoPeak.h:16
tuple node
Definition: Node.py:50
void sparsePeakTreeToStorable(const fftjet::SparseClusteringTree< fftjet::Peak, long > &in, bool writeOutScaleInfo, reco::PattRecoTree< Real, reco::PattRecoPeak< Real > > *out)
void densePeakTreeToStorable(const fftjet::AbsClusteringTree< fftjet::Peak, long > &in, bool writeOutScaleInfo, reco::PattRecoTree< Real, reco::PattRecoPeak< Real > > *out)
tuple out
Definition: dbtoconf.py:99
double Real
void sparsePeakTreeFromStorable(const reco::PattRecoTree< Real, reco::PattRecoPeak< Real > > &in, const std::vector< double > *scaleSetIfNotAdaptive, double completeEventScale, fftjet::SparseClusteringTree< fftjet::Peak, long > *out)
void densePeakTreeFromStorable(const reco::PattRecoTree< Real, reco::PattRecoPeak< Real > > &in, const std::vector< double > *scaleSetIfNotAdaptive, double completeEventScale, fftjet::AbsClusteringTree< fftjet::Peak, long > *out)
edm::TrieNode< PDet > Node