CMS 3D CMS Logo

CMSSW_4_4_3_patch1/src/RecoJets/FFTJetAlgorithms/interface/clusteringTreeConverters.h

Go to the documentation of this file.
00001 //
00002 // Code converting data between FFTJet clustering tree objects and
00003 // event-storable entities. The data can be stored in either single
00004 // or double precision. If the data is stored in double precision,
00005 // the clustering tree can be restored in its original form.
00006 //
00007 // Written by:
00008 //
00009 // Terence Libeiro
00010 // Igor Volobouev
00011 //
00012 // June 2010
00013 
00014 #ifndef RecoJets_FFTJetAlgorithms_clusteringTreeConverters_h
00015 #define RecoJets_FFTJetAlgorithms_clusteringTreeConverters_h
00016 
00017 #include <cfloat>
00018 
00019 #include "fftjet/Peak.hh"
00020 #include "fftjet/AbsClusteringTree.hh"
00021 #include "fftjet/SparseClusteringTree.hh"
00022 
00023 #include "FWCore/Utilities/interface/Exception.h"
00024 #include "DataFormats/JetReco/interface/PattRecoTree.h"
00025 #include "DataFormats/JetReco/interface/PattRecoPeak.h"
00026 
00027 namespace fftjetcms {
00028     // The function below makes PattRecoTree storable in the event
00029     // record out of a sparse clustering tree of fftjet::Peak objects
00030     template<class Real>
00031     void sparsePeakTreeToStorable(
00032         const fftjet::SparseClusteringTree<fftjet::Peak,long>& in,
00033         bool writeOutScaleInfo,
00034         reco::PattRecoTree<Real,reco::PattRecoPeak<Real> >* out);
00035 
00036     // The function below restores sparse clustering tree of fftjet::Peak
00037     // objects using the PattRecoTree. The "completeEventScale" parameter
00038     // should be set to the appropriate scale in case the complete event
00039     // was written out and its scale not stored in the event. If the complete
00040     // event was not written out, the scale must be set to 0.
00041     template<class Real>
00042     void sparsePeakTreeFromStorable(
00043         const reco::PattRecoTree<Real,reco::PattRecoPeak<Real> >& in,
00044         const std::vector<double>* scaleSetIfNotAdaptive,
00045         double completeEventScale,
00046         fftjet::SparseClusteringTree<fftjet::Peak,long>* out);
00047 
00048     // The function below makes PattRecoTree storable in the event
00049     // record out of a dense clustering tree of fftjet::Peak objects
00050     template<class Real>
00051     void densePeakTreeToStorable(
00052         const fftjet::AbsClusteringTree<fftjet::Peak,long>& in,
00053         bool writeOutScaleInfo,
00054         reco::PattRecoTree<Real,reco::PattRecoPeak<Real> >* out);
00055 
00056     // The function below restores dense clustering tree of fftjet::Peak
00057     // objects using the PattRecoTree. The "completeEventScale" parameter
00058     // should be set to the appropriate scale in case the complete event
00059     // was written out and its scale not stored in the event. If the complete
00060     // event was not written out, the scale must be set to 0.
00061     template<class Real>
00062     void densePeakTreeFromStorable(
00063         const reco::PattRecoTree<Real,reco::PattRecoPeak<Real> >& in,
00064         const std::vector<double>* scaleSetIfNotAdaptive,
00065         double completeEventScale,
00066         fftjet::AbsClusteringTree<fftjet::Peak,long>* out);
00067 }
00068 
00070 //
00071 //  Implementation follows
00072 //
00074 
00075 namespace fftjetcms {
00076     // The function below makes PattRecoTree storable in the event
00077     // record out of a sparse clustering tree of fftjet::Peak objects
00078     template<class Real>
00079     void sparsePeakTreeToStorable(
00080         const fftjet::SparseClusteringTree<fftjet::Peak,long>& sparseTree,
00081         const bool writeOutScaleInfo,
00082         reco::PattRecoTree<Real,reco::PattRecoPeak<Real> >* tree)
00083     {
00084         typedef fftjet::SparseClusteringTree<fftjet::Peak,long> SparseTree;
00085         typedef reco::PattRecoPeak<Real> StoredPeak;
00086         typedef reco::PattRecoNode<StoredPeak> StoredNode;
00087         typedef reco::PattRecoTree<Real,StoredPeak> StoredTree;
00088 
00089         assert(tree);
00090 
00091         tree->clear();
00092         tree->setSparse(true);
00093 
00094         const unsigned nNodes = sparseTree.size();
00095         double hessian[3] = {0., 0., 0.};
00096 
00097         // Do not write out the meaningless top node
00098         tree->reserveNodes(nNodes - 1);
00099         for (unsigned i=1; i<nNodes; ++i)
00100         {
00101             const SparseTree::Node& node(sparseTree.getNode(i));
00102             const fftjet::Peak& peak(node.getCluster());
00103             peak.hessian(hessian);
00104             StoredNode sn(StoredPeak(peak.eta(),
00105                                      peak.phi(),
00106                                      peak.magnitude(),
00107                                      hessian,
00108                                      peak.driftSpeed(),
00109                                      peak.magSpeed(),
00110                                      peak.lifetime(),
00111                                      peak.scale(),
00112                                      peak.nearestNeighborDistance(),
00113                                      peak.clusterRadius(),
00114                                      peak.clusterSeparation()),
00115                           node.originalLevel(),
00116                           node.mask(),
00117                           node.parent());
00118             tree->addNode(sn);
00119         }
00120 
00121         // Do we want to write out the scales? We will use the following
00122         // convention: if the tree is using an adaptive algorithm, the scales
00123         // will be written out. If not, they are not going to change from
00124         // event to event. In this case the scales would waste disk space
00125         // for no particularly good reason, so they will not be written out.
00126         if (writeOutScaleInfo)
00127         {
00128             // Do not write out the meaningless top-level scale
00129             const unsigned nScales = sparseTree.nLevels();
00130             tree->reserveScales(nScales - 1);
00131             for (unsigned i=1; i<nScales; ++i)
00132                 tree->addScale(sparseTree.getScale(i));
00133         }
00134     }
00135 
00136     // The function below restores sparse clustering tree of fftjet::Peak
00137     // objects using the PattRecoTree
00138     template<class Real>
00139     void sparsePeakTreeFromStorable(
00140         const reco::PattRecoTree<Real,reco::PattRecoPeak<Real> >& in,
00141         const std::vector<double>* scaleSetIfNotAdaptive,
00142         const double completeEventScale,
00143         fftjet::SparseClusteringTree<fftjet::Peak,long>* out)
00144     {
00145         typedef fftjet::SparseClusteringTree<fftjet::Peak,long> SparseTree;
00146         typedef reco::PattRecoPeak<Real> StoredPeak;
00147         typedef reco::PattRecoNode<StoredPeak> StoredNode;
00148         typedef reco::PattRecoTree<Real,StoredPeak> StoredTree;
00149 
00150         if (!in.isSparse())
00151             throw cms::Exception("FFTJetBadConfig")
00152                 << "can't restore sparse clustering tree"
00153                 << " from densely stored record" << std::endl;
00154 
00155         assert(out);
00156         out->clear();
00157 
00158         const std::vector<StoredNode>& nodes(in.getNodes());
00159         const unsigned n = nodes.size();
00160         out->reserveNodes(n + 1U);
00161 
00162         double hessian[3] = {0., 0., 0.};
00163 
00164         for (unsigned i=0; i<n; ++i)
00165         {
00166             const StoredNode& snode(nodes[i]);
00167             const StoredPeak& p(snode.getCluster());
00168             p.hessian(hessian);
00169             const SparseTree::Node node(
00170                 fftjet::Peak(p.eta(), p.phi(), p.magnitude(),
00171                              hessian, p.driftSpeed(),
00172                              p.magSpeed(), p.lifetime(),
00173                              p.scale(), p.nearestNeighborDistance(),
00174                              1.0, 0.0, 0.0,
00175                              p.clusterRadius(), p.clusterSeparation()),
00176                 snode.originalLevel(),
00177                 snode.mask());
00178             out->addNode(node, snode.parent());
00179         }
00180 
00181         const std::vector<Real>& storedScales(in.getScales());
00182         if (!storedScales.empty())
00183         {
00184             const unsigned nsc = storedScales.size();
00185             out->reserveScales(nsc + 1U);
00186             out->addScale(DBL_MAX);
00187             const Real* scales = &storedScales[0];
00188             for (unsigned i=0; i<nsc; ++i)
00189                 out->addScale(scales[i]);
00190         }
00191         else if (scaleSetIfNotAdaptive && !scaleSetIfNotAdaptive->empty())
00192         {
00193             const unsigned nsc = scaleSetIfNotAdaptive->size();
00194             // There may be the "complete event" scale added at the end.
00195             // Reserve a sufficient number of scales to take this into
00196             // account.
00197             if (completeEventScale)
00198                 out->reserveScales(nsc + 2U);
00199             else
00200                 out->reserveScales(nsc + 1U);
00201             out->addScale(DBL_MAX);
00202             const double* scales = &(*scaleSetIfNotAdaptive)[0];
00203             for (unsigned i=0; i<nsc; ++i)
00204                 out->addScale(scales[i]);
00205             if (completeEventScale)
00206                 out->addScale(completeEventScale);
00207         }
00208         else
00209         {
00210             throw cms::Exception("FFTJetBadConfig")
00211                 << "can't restore sparse clustering tree scales"
00212                 << std::endl;
00213         }
00214     }
00215 
00216     // The function below makes PattRecoTree storable in the event
00217     // record out of a dense clustering tree of fftjet::Peak objects
00218     template<class Real>
00219     void densePeakTreeToStorable(
00220         const fftjet::AbsClusteringTree<fftjet::Peak,long>& in,
00221         bool writeOutScaleInfo,
00222         reco::PattRecoTree<Real,reco::PattRecoPeak<Real> >* out)
00223     {
00224         typedef fftjet::AbsClusteringTree<fftjet::Peak,long> DenseTree;
00225         typedef reco::PattRecoPeak<Real> StoredPeak;
00226         typedef reco::PattRecoNode<StoredPeak> StoredNode;
00227         typedef reco::PattRecoTree<Real,StoredPeak> StoredTree;
00228         
00229         
00230         assert(out);
00231         out->clear();
00232         out->setSparse(false);
00233         
00234         
00235         const unsigned nLevels = in.nLevels();
00236         double hessian[3] = {0., 0., 0.};
00237 
00238         
00239         // Do not write out the meaningless top node
00240         out->reserveNodes(in.nClusters() - 1);
00241         
00242         
00243         for (unsigned i=1; i<nLevels; ++i)
00244         {
00245             
00246             const unsigned int nclus = in.nClusters(i);
00247             DenseTree::NodeId id(i,0);
00248             for (;id.second<nclus; ++id.second) 
00249             {
00250            
00251            
00252                 const fftjet::Peak& peak(in.getCluster(id));
00253                 peak.hessian(hessian);
00254                 StoredNode sn(StoredPeak(peak.eta(),
00255                                          peak.phi(),
00256                                          peak.magnitude(),
00257                                          hessian,
00258                                          peak.driftSpeed(),
00259                                          peak.magSpeed(),
00260                                          peak.lifetime(),
00261                                          peak.scale(),
00262                                          peak.nearestNeighborDistance(),
00263                                          peak.clusterRadius(),
00264                                          peak.clusterSeparation()),
00265                               i,
00266                               0,
00267                               0);
00268                 
00269                 out->addNode(sn);
00270                 
00271             }
00272             
00273         }
00274         
00275     
00276         // Do we want to write out the scales? We will use the following
00277         // convention: if the tree is using an adaptive algorithm, the scales
00278         // will be written out. If not, they are not going to change from
00279         // event to event. In this case the scales would waste disk space
00280         // for no particularly good reason, so they will not be written out.
00281         if (writeOutScaleInfo)
00282         {
00283             // Do not write out the meaningless top-level scale
00284             const unsigned nScales = in.nLevels();
00285             
00286             out->reserveScales(nScales - 1);
00287             
00288             for (unsigned i=1; i<nScales; ++i)
00289                 out->addScale(in.getScale(i));
00290             
00291         }
00292     }
00293 
00294     // The function below restores dense clustering tree of fftjet::Peak
00295     // objects using the PattRecoTree
00296     template<class Real>
00297     void densePeakTreeFromStorable(
00298         const reco::PattRecoTree<Real,reco::PattRecoPeak<Real> >& in,
00299         const std::vector<double>* scaleSetIfNotAdaptive,
00300         const double completeEventScale,
00301         fftjet::AbsClusteringTree<fftjet::Peak,long>* out)
00302     {
00303         typedef fftjet::AbsClusteringTree<fftjet::Peak,long> DenseTree;
00304         typedef reco::PattRecoPeak<Real> StoredPeak;
00305         typedef reco::PattRecoNode<StoredPeak> StoredNode;
00306         typedef reco::PattRecoTree<Real,StoredPeak> StoredTree;
00307 
00308         if (in.isSparse())
00309             throw cms::Exception("FFTJetBadConfig")
00310                 << "can't restore dense clustering tree"
00311                 << " from sparsely stored record" << std::endl;
00312 
00313         assert(out);
00314         out->clear();
00315 
00316         const std::vector<StoredNode>& nodes(in.getNodes());
00317         const unsigned n = nodes.size();
00318         double hessian[3] = {0., 0., 0.};
00319                 
00320         const std::vector<Real>& scales (in.getScales());
00321         unsigned int scnum     = 0;
00322         std::vector<fftjet::Peak> clusters;
00323         const unsigned scsize1 = scales.size();
00324         unsigned scsize2 = scaleSetIfNotAdaptive ? scaleSetIfNotAdaptive->size() : 0;
00325         if (scsize2 && completeEventScale) ++scsize2;
00326         const unsigned scsize  = (scsize1==0?scsize2:scsize1);
00327 
00328         if (scsize == 0)  
00329             throw cms::Exception("FFTJetBadConfig")
00330                 << " No scales passed to the function densePeakTreeFromStorable()"
00331                 << std::endl;
00332 
00333         // to check whether the largest level equals the size of scale vector
00334         const double* sc_not_ad = scsize2 ? &(*scaleSetIfNotAdaptive)[0] : 0;
00335 
00336         unsigned templevel = n ? nodes[0].originalLevel() : 1;
00337         for (unsigned i=0; i<n; ++i)
00338         {
00339             const StoredNode& snode(nodes[i]);
00340             const StoredPeak& p(snode.getCluster());
00341             p.hessian(hessian);
00342         
00343             const unsigned levelNumber = snode.originalLevel();
00344 
00345             if (templevel != levelNumber) 
00346             {
00347                 if (scnum >= scsize)
00348                     throw cms::Exception("FFTJetBadConfig")
00349                         << "bad scales, please check the scales"
00350                         << std::endl;
00351                 const double scale = ( (scsize1==0) ? sc_not_ad[scnum] : scales[scnum] );
00352                 out->insert(scale, clusters, 0L);
00353                 clusters.clear();
00354                 templevel = levelNumber;
00355                 ++scnum;
00356             }
00357 
00358             fftjet::Peak apeak(p.eta(), p.phi(), p.magnitude(),
00359                                hessian, p.driftSpeed(),
00360                                p.magSpeed(), p.lifetime(),
00361                                p.scale(), p.nearestNeighborDistance(),
00362                                1.0, 0.0, 0.0,
00363                                p.clusterRadius(), p.clusterSeparation());
00364             clusters.push_back(apeak);
00365                
00366             if (i==(n-1) && levelNumber!=scsize) 
00367                 throw cms::Exception("FFTJetBadConfig")
00368                     << "bad scales, please check the scales"
00369                     << std::endl;
00370         }
00371 
00372         const double scale = scsize1 ? scales[scnum] : completeEventScale ? 
00373             completeEventScale : sc_not_ad[scnum];
00374         out->insert(scale, clusters, 0L);
00375     }
00376 }
00377 
00378 #endif // RecoJets_FFTJetAlgorithms_clusteringTreeConverters_h