14 #ifndef RecoJets_FFTJetAlgorithms_clusteringTreeConverters_h
15 #define RecoJets_FFTJetAlgorithms_clusteringTreeConverters_h
19 #include "fftjet/Peak.hh"
20 #include "fftjet/AbsClusteringTree.hh"
21 #include "fftjet/SparseClusteringTree.hh"
32 bool writeOutScaleInfo,
42 const std::vector<double>* scaleSetIfNotAdaptive,
44 fftjet::SparseClusteringTree<fftjet::Peak, long>*
out);
50 bool writeOutScaleInfo,
60 const std::vector<double>* scaleSetIfNotAdaptive,
62 fftjet::AbsClusteringTree<fftjet::Peak, long>*
out);
76 const bool writeOutScaleInfo,
78 typedef fftjet::SparseClusteringTree<fftjet::Peak, long> SparseTree;
86 tree->setSparse(
true);
88 const unsigned nNodes = sparseTree.size();
89 double hessian[3] = {0., 0., 0.};
92 tree->reserveNodes(nNodes - 1);
93 for (
unsigned i = 1;
i < nNodes; ++
i) {
95 const fftjet::Peak& peak(node.getCluster());
96 peak.hessian(hessian);
97 StoredNode sn(StoredPeak(peak.eta(),
105 peak.nearestNeighborDistance(),
106 peak.clusterRadius(),
107 peak.clusterSeparation(),
110 node.originalLevel(),
121 if (writeOutScaleInfo) {
123 const unsigned nScales = sparseTree.nLevels();
126 tree->addScale(sparseTree.getScale(
i));
132 template <
class Real>
134 const std::vector<double>* scaleSetIfNotAdaptive,
136 fftjet::SparseClusteringTree<fftjet::Peak, long>*
out) {
137 typedef fftjet::SparseClusteringTree<fftjet::Peak, long> SparseTree;
143 throw cms::Exception(
"FFTJetBadConfig") <<
"can't restore sparse clustering tree"
144 <<
" from densely stored record" << std::endl;
149 const std::vector<StoredNode>&
nodes(
in.getNodes());
150 const unsigned n =
nodes.size();
151 out->reserveNodes(
n + 1
U);
153 double hessian[3] = {0., 0., 0.};
155 for (
unsigned i = 0;
i <
n; ++
i) {
156 const StoredNode& snode(
nodes[
i]);
157 const StoredPeak&
p(snode.getCluster());
159 fftjet::Peak peak(
p.eta(),
167 p.nearestNeighborDistance(),
172 p.clusterSeparation());
173 peak.setSplitTime(
p.splitTime());
174 peak.setMergeTime(
p.mergeTime());
176 out->addNode(node, snode.parent());
179 const std::vector<Real>& storedScales(
in.getScales());
180 if (!storedScales.empty()) {
181 const unsigned nsc = storedScales.size();
182 out->reserveScales(nsc + 1
U);
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();
193 out->reserveScales(nsc + 2
U);
195 out->reserveScales(nsc + 1
U);
196 out->addScale(DBL_MAX);
197 const double* scales = &(*scaleSetIfNotAdaptive)[0];
198 for (
unsigned i = 0;
i < nsc; ++
i)
199 out->addScale(scales[
i]);
203 throw cms::Exception(
"FFTJetBadConfig") <<
"can't restore sparse clustering tree scales" << std::endl;
209 template <
class Real>
211 bool writeOutScaleInfo,
213 typedef fftjet::AbsClusteringTree<fftjet::Peak, long> DenseTree;
220 out->setSparse(
false);
222 const unsigned nLevels =
in.nLevels();
223 double hessian[3] = {0., 0., 0.};
226 out->reserveNodes(
in.nClusters() - 1);
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(),
242 peak.nearestNeighborDistance(),
243 peak.clusterRadius(),
244 peak.clusterSeparation(),
259 if (writeOutScaleInfo) {
264 out->addScale(
in.getScale(
i));
270 template <
class Real>
272 const std::vector<double>* scaleSetIfNotAdaptive,
274 fftjet::AbsClusteringTree<fftjet::Peak, long>*
out) {
275 typedef fftjet::AbsClusteringTree<fftjet::Peak, long> DenseTree;
281 throw cms::Exception(
"FFTJetBadConfig") <<
"can't restore dense clustering tree"
282 <<
" from sparsely stored record" << std::endl;
287 const std::vector<StoredNode>&
nodes(
in.getNodes());
288 const unsigned n =
nodes.size();
289 double hessian[3] = {0., 0., 0.};
291 const std::vector<Real>& scales(
in.getScales());
292 unsigned int scnum = 0;
294 const unsigned scsize1 = scales.size();
295 unsigned scsize2 = scaleSetIfNotAdaptive ? scaleSetIfNotAdaptive->size() : 0;
298 const unsigned scsize = (scsize1 == 0 ? scsize2 : scsize1);
302 <<
" No scales passed to the function densePeakTreeFromStorable()" << std::endl;
305 const double* sc_not_ad = scsize2 ? &(*scaleSetIfNotAdaptive)[0] :
nullptr;
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());
313 const unsigned levelNumber = snode.originalLevel();
315 if (templevel != levelNumber) {
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]);
321 templevel = levelNumber;
325 fftjet::Peak apeak(
p.eta(),
333 p.nearestNeighborDistance(),
338 p.clusterSeparation());
339 apeak.setSplitTime(
p.splitTime());
340 apeak.setMergeTime(
p.mergeTime());
343 if (
i == (
n - 1) && levelNumber != scsize)
344 throw cms::Exception(
"FFTJetBadConfig") <<
"bad scales, please check the scales" << std::endl;
352 #endif // RecoJets_FFTJetAlgorithms_clusteringTreeConverters_h