CMS 3D CMS Logo

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>
31  void sparsePeakTreeToStorable(const fftjet::SparseClusteringTree<fftjet::Peak, long>& in,
32  bool writeOutScaleInfo,
34 
35  // The function below restores sparse clustering tree of fftjet::Peak
36  // objects using the PattRecoTree. The "completeEventScale" parameter
37  // should be set to the appropriate scale in case the complete event
38  // was written out and its scale not stored in the event. If the complete
39  // event was not written out, the scale must be set to 0.
40  template <class Real>
42  const std::vector<double>* scaleSetIfNotAdaptive,
43  double completeEventScale,
44  fftjet::SparseClusteringTree<fftjet::Peak, long>* out);
45 
46  // The function below makes PattRecoTree storable in the event
47  // record out of a dense clustering tree of fftjet::Peak objects
48  template <class Real>
49  void densePeakTreeToStorable(const fftjet::AbsClusteringTree<fftjet::Peak, long>& in,
50  bool writeOutScaleInfo,
52 
53  // The function below restores dense clustering tree of fftjet::Peak
54  // objects using the PattRecoTree. The "completeEventScale" parameter
55  // should be set to the appropriate scale in case the complete event
56  // was written out and its scale not stored in the event. If the complete
57  // event was not written out, the scale must be set to 0.
58  template <class Real>
60  const std::vector<double>* scaleSetIfNotAdaptive,
61  double completeEventScale,
62  fftjet::AbsClusteringTree<fftjet::Peak, long>* out);
63 } // namespace fftjetcms
64 
66 //
67 // Implementation follows
68 //
70 
71 namespace fftjetcms {
72  // The function below makes PattRecoTree storable in the event
73  // record out of a sparse clustering tree of fftjet::Peak objects
74  template <class Real>
75  void sparsePeakTreeToStorable(const fftjet::SparseClusteringTree<fftjet::Peak, long>& sparseTree,
76  const bool writeOutScaleInfo,
78  typedef fftjet::SparseClusteringTree<fftjet::Peak, long> SparseTree;
79  typedef reco::PattRecoPeak<Real> StoredPeak;
80  typedef reco::PattRecoNode<StoredPeak> StoredNode;
81  typedef reco::PattRecoTree<Real, StoredPeak> StoredTree;
82 
83  assert(tree);
84 
85  tree->clear();
86  tree->setSparse(true);
87 
88  const unsigned nNodes = sparseTree.size();
89  double hessian[3] = {0., 0., 0.};
90 
91  // Do not write out the meaningless top node
92  tree->reserveNodes(nNodes - 1);
93  for (unsigned i = 1; i < nNodes; ++i) {
94  const SparseTree::Node& node(sparseTree.getNode(i));
95  const fftjet::Peak& peak(node.getCluster());
96  peak.hessian(hessian);
97  StoredNode sn(StoredPeak(peak.eta(),
98  peak.phi(),
99  peak.magnitude(),
100  hessian,
101  peak.driftSpeed(),
102  peak.magSpeed(),
103  peak.lifetime(),
104  peak.scale(),
105  peak.nearestNeighborDistance(),
106  peak.clusterRadius(),
107  peak.clusterSeparation(),
108  peak.splitTime(),
109  peak.mergeTime()),
110  node.originalLevel(),
111  node.mask(),
112  node.parent());
113  tree->addNode(sn);
114  }
115 
116  // Do we want to write out the scales? We will use the following
117  // convention: if the tree is using an adaptive algorithm, the scales
118  // will be written out. If not, they are not going to change from
119  // event to event. In this case the scales would waste disk space
120  // for no particularly good reason, so they will not be written out.
121  if (writeOutScaleInfo) {
122  // Do not write out the meaningless top-level scale
123  const unsigned nScales = sparseTree.nLevels();
124  tree->reserveScales(nScales - 1);
125  for (unsigned i = 1; i < nScales; ++i)
126  tree->addScale(sparseTree.getScale(i));
127  }
128  }
129 
130  // The function below restores sparse clustering tree of fftjet::Peak
131  // objects using the PattRecoTree
132  template <class Real>
134  const std::vector<double>* scaleSetIfNotAdaptive,
135  const double completeEventScale,
136  fftjet::SparseClusteringTree<fftjet::Peak, long>* out) {
137  typedef fftjet::SparseClusteringTree<fftjet::Peak, long> SparseTree;
138  typedef reco::PattRecoPeak<Real> StoredPeak;
139  typedef reco::PattRecoNode<StoredPeak> StoredNode;
140  typedef reco::PattRecoTree<Real, StoredPeak> StoredTree;
141 
142  if (!in.isSparse())
143  throw cms::Exception("FFTJetBadConfig") << "can't restore sparse clustering tree"
144  << " from densely stored record" << std::endl;
145 
146  assert(out);
147  out->clear();
148 
149  const std::vector<StoredNode>& nodes(in.getNodes());
150  const unsigned n = nodes.size();
151  out->reserveNodes(n + 1U);
152 
153  double hessian[3] = {0., 0., 0.};
154 
155  for (unsigned i = 0; i < n; ++i) {
156  const StoredNode& snode(nodes[i]);
157  const StoredPeak& p(snode.getCluster());
158  p.hessian(hessian);
159  fftjet::Peak peak(p.eta(),
160  p.phi(),
161  p.magnitude(),
162  hessian,
163  p.driftSpeed(),
164  p.magSpeed(),
165  p.lifetime(),
166  p.scale(),
167  p.nearestNeighborDistance(),
168  1.0,
169  0.0,
170  0.0,
171  p.clusterRadius(),
172  p.clusterSeparation());
173  peak.setSplitTime(p.splitTime());
174  peak.setMergeTime(p.mergeTime());
175  const SparseTree::Node node(peak, snode.originalLevel(), snode.mask());
176  out->addNode(node, snode.parent());
177  }
178 
179  const std::vector<Real>& storedScales(in.getScales());
180  if (!storedScales.empty()) {
181  const unsigned nsc = storedScales.size();
182  out->reserveScales(nsc + 1U);
183  out->addScale(DBL_MAX);
184  const Real* scales = &storedScales[0];
185  for (unsigned i = 0; i < nsc; ++i)
186  out->addScale(scales[i]);
187  } else if (scaleSetIfNotAdaptive && !scaleSetIfNotAdaptive->empty()) {
188  const unsigned nsc = scaleSetIfNotAdaptive->size();
189  // There may be the "complete event" scale added at the end.
190  // Reserve a sufficient number of scales to take this into
191  // account.
192  if (completeEventScale)
193  out->reserveScales(nsc + 2U);
194  else
195  out->reserveScales(nsc + 1U);
196  out->addScale(DBL_MAX);
197  const double* scales = &(*scaleSetIfNotAdaptive)[0];
198  for (unsigned i = 0; i < nsc; ++i)
199  out->addScale(scales[i]);
200  if (completeEventScale)
201  out->addScale(completeEventScale);
202  } else {
203  throw cms::Exception("FFTJetBadConfig") << "can't restore sparse clustering tree scales" << std::endl;
204  }
205  }
206 
207  // The function below makes PattRecoTree storable in the event
208  // record out of a dense clustering tree of fftjet::Peak objects
209  template <class Real>
210  void densePeakTreeToStorable(const fftjet::AbsClusteringTree<fftjet::Peak, long>& in,
211  bool writeOutScaleInfo,
213  typedef fftjet::AbsClusteringTree<fftjet::Peak, long> DenseTree;
214  typedef reco::PattRecoPeak<Real> StoredPeak;
215  typedef reco::PattRecoNode<StoredPeak> StoredNode;
216  typedef reco::PattRecoTree<Real, StoredPeak> StoredTree;
217 
218  assert(out);
219  out->clear();
220  out->setSparse(false);
221 
222  const unsigned nLevels = in.nLevels();
223  double hessian[3] = {0., 0., 0.};
224 
225  // Do not write out the meaningless top node
226  out->reserveNodes(in.nClusters() - 1);
227 
228  for (unsigned i = 1; i < nLevels; ++i) {
229  const unsigned int nclus = in.nClusters(i);
230  DenseTree::NodeId id(i, 0);
231  for (; id.second < nclus; ++id.second) {
232  const fftjet::Peak& peak(in.getCluster(id));
233  peak.hessian(hessian);
234  StoredNode sn(StoredPeak(peak.eta(),
235  peak.phi(),
236  peak.magnitude(),
237  hessian,
238  peak.driftSpeed(),
239  peak.magSpeed(),
240  peak.lifetime(),
241  peak.scale(),
242  peak.nearestNeighborDistance(),
243  peak.clusterRadius(),
244  peak.clusterSeparation(),
245  peak.splitTime(),
246  peak.mergeTime()),
247  i,
248  0,
249  0);
250  out->addNode(sn);
251  }
252  }
253 
254  // Do we want to write out the scales? We will use the following
255  // convention: if the tree is using an adaptive algorithm, the scales
256  // will be written out. If not, they are not going to change from
257  // event to event. In this case the scales would waste disk space
258  // for no particularly good reason, so they will not be written out.
259  if (writeOutScaleInfo) {
260  // Do not write out the meaningless top-level scale
261  const unsigned nScales = in.nLevels();
262  out->reserveScales(nScales - 1);
263  for (unsigned i = 1; i < nScales; ++i)
264  out->addScale(in.getScale(i));
265  }
266  }
267 
268  // The function below restores dense clustering tree of fftjet::Peak
269  // objects using the PattRecoTree
270  template <class Real>
272  const std::vector<double>* scaleSetIfNotAdaptive,
273  const double completeEventScale,
274  fftjet::AbsClusteringTree<fftjet::Peak, long>* out) {
275  typedef fftjet::AbsClusteringTree<fftjet::Peak, long> DenseTree;
276  typedef reco::PattRecoPeak<Real> StoredPeak;
277  typedef reco::PattRecoNode<StoredPeak> StoredNode;
278  typedef reco::PattRecoTree<Real, StoredPeak> StoredTree;
279 
280  if (in.isSparse())
281  throw cms::Exception("FFTJetBadConfig") << "can't restore dense clustering tree"
282  << " from sparsely stored record" << std::endl;
283 
284  assert(out);
285  out->clear();
286 
287  const std::vector<StoredNode>& nodes(in.getNodes());
288  const unsigned n = nodes.size();
289  double hessian[3] = {0., 0., 0.};
290 
291  const std::vector<Real>& scales(in.getScales());
292  unsigned int scnum = 0;
293  std::vector<fftjet::Peak> clusters;
294  const unsigned scsize1 = scales.size();
295  unsigned scsize2 = scaleSetIfNotAdaptive ? scaleSetIfNotAdaptive->size() : 0;
296  if (scsize2 && completeEventScale)
297  ++scsize2;
298  const unsigned scsize = (scsize1 == 0 ? scsize2 : scsize1);
299 
300  if (scsize == 0)
301  throw cms::Exception("FFTJetBadConfig")
302  << " No scales passed to the function densePeakTreeFromStorable()" << std::endl;
303 
304  // to check whether the largest level equals the size of scale vector
305  const double* sc_not_ad = scsize2 ? &(*scaleSetIfNotAdaptive)[0] : nullptr;
306 
307  unsigned templevel = n ? nodes[0].originalLevel() : 1;
308  for (unsigned i = 0; i < n; ++i) {
309  const StoredNode& snode(nodes[i]);
310  const StoredPeak& p(snode.getCluster());
311  p.hessian(hessian);
312 
313  const unsigned levelNumber = snode.originalLevel();
314 
315  if (templevel != levelNumber) {
316  if (scnum >= scsize)
317  throw cms::Exception("FFTJetBadConfig") << "bad scales, please check the scales" << std::endl;
318  const double scale = ((scsize1 == 0) ? sc_not_ad[scnum] : scales[scnum]);
319  out->insert(scale, clusters, 0L);
320  clusters.clear();
321  templevel = levelNumber;
322  ++scnum;
323  }
324 
325  fftjet::Peak apeak(p.eta(),
326  p.phi(),
327  p.magnitude(),
328  hessian,
329  p.driftSpeed(),
330  p.magSpeed(),
331  p.lifetime(),
332  p.scale(),
333  p.nearestNeighborDistance(),
334  1.0,
335  0.0,
336  0.0,
337  p.clusterRadius(),
338  p.clusterSeparation());
339  apeak.setSplitTime(p.splitTime());
340  apeak.setMergeTime(p.mergeTime());
341  clusters.push_back(apeak);
342 
343  if (i == (n - 1) && levelNumber != scsize)
344  throw cms::Exception("FFTJetBadConfig") << "bad scales, please check the scales" << std::endl;
345  }
346 
347  const double scale = scsize1 ? scales[scnum] : completeEventScale ? completeEventScale : sc_not_ad[scnum];
348  out->insert(scale, clusters, 0L);
349  }
350 } // namespace fftjetcms
351 
352 #endif // RecoJets_FFTJetAlgorithms_clusteringTreeConverters_h
fftjetcms::densePeakTreeToStorable
void densePeakTreeToStorable(const fftjet::AbsClusteringTree< fftjet::Peak, long > &in, bool writeOutScaleInfo, reco::PattRecoTree< Real, reco::PattRecoPeak< Real > > *out)
Definition: clusteringTreeConverters.h:210
dttmaxenums::L
Definition: DTTMax.h:29
mps_fire.i
i
Definition: mps_fire.py:428
dqmiodumpmetadata.n
n
Definition: dqmiodumpmetadata.py:28
fftjetcms
Definition: AbsPileupCalculator.h:15
L1EGammaCrystalsEmulatorProducer_cfi.scale
scale
Definition: L1EGammaCrystalsEmulatorProducer_cfi.py:10
fftjetcommon_cfi.nScales
nScales
Definition: fftjetcommon_cfi.py:111
cms::Node
TGeoNode Node
Definition: DDFilteredView.h:55
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
tree
Definition: tree.py:1
cms::cuda::assert
assert(be >=bs)
fftjetdijetfilter_cfi.completeEventScale
completeEventScale
Definition: fftjetdijetfilter_cfi.py:17
fftjetcms::densePeakTreeFromStorable
void densePeakTreeFromStorable(const reco::PattRecoTree< Real, reco::PattRecoPeak< Real > > &in, const std::vector< double > *scaleSetIfNotAdaptive, double completeEventScale, fftjet::AbsClusteringTree< fftjet::Peak, long > *out)
Definition: clusteringTreeConverters.h:271
fftjetcms::sparsePeakTreeToStorable
void sparsePeakTreeToStorable(const fftjet::SparseClusteringTree< fftjet::Peak, long > &in, bool writeOutScaleInfo, reco::PattRecoTree< Real, reco::PattRecoPeak< Real > > *out)
Definition: clusteringTreeConverters.h:75
class-composition.nodes
nodes
Definition: class-composition.py:75
fftjetcms::sparsePeakTreeFromStorable
void sparsePeakTreeFromStorable(const reco::PattRecoTree< Real, reco::PattRecoPeak< Real > > &in, const std::vector< double > *scaleSetIfNotAdaptive, double completeEventScale, fftjet::SparseClusteringTree< fftjet::Peak, long > *out)
Definition: clusteringTreeConverters.h:133
mitigatedMETSequence_cff.U
U
Definition: mitigatedMETSequence_cff.py:36
fftjetcms::Real
double Real
Definition: fftjetTypedefs.h:21
bsc_activity_cfg.clusters
clusters
Definition: bsc_activity_cfg.py:36
reco::PattRecoTree
Class for storing FFTJet sparse clustering trees.
Definition: PattRecoTree.h:20
reco::PattRecoPeak
Preclusters from FFTJet pattern recognition stage.
Definition: PattRecoPeak.h:16
recoMuon::in
Definition: RecoMuonEnumerators.h:6
reco::PattRecoNode
Tree nodes for storing FFTJet preclusters.
Definition: PattRecoNode.h:16
PattRecoTree.h
triggerObjects_cff.id
id
Definition: triggerObjects_cff.py:29
Exception
Definition: hltDiff.cc:245
Exception.h
MillePedeFileConverter_cfg.out
out
Definition: MillePedeFileConverter_cfg.py:31
cms::Exception
Definition: Exception.h:70
PattRecoPeak.h