CMS 3D CMS Logo

FFTJetTreeDump.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: FFTJetTreeDump
4 // Class: FFTJetTreeDump
5 //
13 //
14 // Original Author: Igor Volobouev
15 // Created: Sun Jun 20 14:32:36 CDT 2010
16 //
17 //
18 
19 #include <iostream>
20 #include <sstream>
21 #include <fstream>
22 #include <functional>
23 
24 // FFTJet headers
25 #include "fftjet/ProximityClusteringTree.hh"
26 #include "fftjet/OpenDXPeakTree.hh"
27 
28 // user include files
31 
34 
36 
37 // parameter parser header
39 
40 // functions which manipulate storable trees
42 
44 
45 using namespace fftjetcms;
46 
47 //
48 // class declaration
49 //
51 public:
52  explicit FFTJetTreeDump(const edm::ParameterSet&);
53  ~FFTJetTreeDump() override;
54 
55 private:
56  // Useful local typedefs
57  typedef fftjet::ProximityClusteringTree<fftjet::Peak, long> ClusteringTree;
58  typedef fftjet::SparseClusteringTree<fftjet::Peak, long> SparseTree;
59  typedef fftjet::OpenDXPeakTree<long, fftjet::AbsClusteringTree> DXFormatter;
60  typedef fftjet::OpenDXPeakTree<long, fftjet::SparseClusteringTree> SparseFormatter;
61  typedef fftjet::Functor1<double, fftjet::Peak> PeakProperty;
63 
64  FFTJetTreeDump() = delete;
65  FFTJetTreeDump(const FFTJetTreeDump&) = delete;
66  FFTJetTreeDump& operator=(const FFTJetTreeDump&) = delete;
67 
68  void beginJob() override;
69  void analyze(const edm::Event&, const edm::EventSetup&) override;
70  void endJob() override;
71 
72  template <class Real>
73  void processTreeData(const edm::Event&, std::ofstream&);
74 
75  template <class Ptr>
76  void checkConfig(const Ptr& ptr, const char* message) {
77  if (ptr.get() == nullptr)
78  throw cms::Exception("FFTJetBadConfig") << message << std::endl;
79  }
80 
81  // ----------member data ---------------------------
82  // The complete clustering tree
83  ClusteringTree* clusteringTree;
84 
87 
89  const double etaMax;
91  const bool insertCompleteEvent;
92  const double completeEventScale;
93 
94  // Distance calculator for the clustering tree
95  std::unique_ptr<fftjet::AbsDistanceCalculator<fftjet::Peak> > distanceCalc;
96 
97  // Scales used
98  std::unique_ptr<std::vector<double> > iniScales;
99 
100  // The sparse clustering tree
101  SparseTree sparseTree;
102 
103  // Functors which define OpenDX glyph size and color
104  std::unique_ptr<PeakProperty> glyphSize;
105  std::unique_ptr<PeakProperty> glyphColor;
106 
107  // OpenDX formatters
108  std::unique_ptr<DXFormatter> denseFormatter;
109  std::unique_ptr<SparseFormatter> sparseFormatter;
110 
111  unsigned counter;
112 };
113 
114 //
115 // constructors and destructor
116 //
118  : clusteringTree(nullptr),
119  treeLabel(ps.getParameter<edm::InputTag>("treeLabel")),
120  outputPrefix(ps.getParameter<std::string>("outputPrefix")),
121  etaMax(ps.getParameter<double>("etaMax")),
122  storeInSinglePrecision(true),
123  insertCompleteEvent(ps.getParameter<bool>("insertCompleteEvent")),
124  completeEventScale(ps.getParameter<double>("completeEventScale")),
125  counter(0) {
126  if (etaMax < 0.0)
127  throw cms::Exception("FFTJetBadConfig") << "etaMax can not be negative" << std::endl;
128 
129  // Build the set of pattern recognition scales
130  const edm::ParameterSet& InitialScales(ps.getParameter<edm::ParameterSet>("InitialScales"));
131  iniScales = fftjet_ScaleSet_parser(InitialScales);
132  checkConfig(iniScales, "invalid set of scales");
133  std::sort(iniScales->begin(), iniScales->end(), std::greater<double>());
134 
135  // Distance calculator for the clustering tree
136  const edm::ParameterSet& TreeDistanceCalculator(ps.getParameter<edm::ParameterSet>("TreeDistanceCalculator"));
137  distanceCalc = fftjet_DistanceCalculator_parser(TreeDistanceCalculator);
138  checkConfig(distanceCalc, "invalid tree distance calculator");
139 
140  // Determine representations for the OpenDX glyph size and color
143  checkConfig(glyphSize, "invalid glyph size parameters");
144 
147  checkConfig(glyphColor, "invalid glyph color parameters");
148 
149  // Build the tree formatters
150  denseFormatter = std::unique_ptr<DXFormatter>(new DXFormatter(glyphSize.get(), glyphColor.get(), etaMax));
151  sparseFormatter = std::unique_ptr<SparseFormatter>(new SparseFormatter(glyphSize.get(), glyphColor.get(), etaMax));
152 
153  // Build the clustering tree
155 
156  treeToken = consumes<StoredTree>(treeLabel);
157 }
158 
160 
161 //
162 // member functions
163 //
164 template <class Real>
166  // Get the event number
167  edm::RunNumber_t const runNum = iEvent.id().run();
168  edm::EventNumber_t const evNum = iEvent.id().event();
169 
170  // Get the input
172  iEvent.getByToken(treeToken, input);
173 
174  const double eventScale = insertCompleteEvent ? completeEventScale : 0.0;
175  if (input->isSparse()) {
176  sparsePeakTreeFromStorable(*input, iniScales.get(), eventScale, &sparseTree);
177  sparseFormatter->setTree(sparseTree, runNum, evNum);
178  file << *sparseFormatter << std::endl;
179  } else {
180  densePeakTreeFromStorable(*input, iniScales.get(), eventScale, clusteringTree);
181  denseFormatter->setTree(*clusteringTree, runNum, evNum);
182  file << *denseFormatter << std::endl;
183  }
184 }
185 
186 // ------------ method called to for each event ------------
188  // Create the file name
189  std::ostringstream filename;
190  filename << outputPrefix << '_' << counter++ << ".dx";
191 
192  // Open the file
193  std::ofstream file(filename.str().c_str());
194  if (!file)
195  throw cms::Exception("FFTJetBadConfig") << "Failed to open file \"" << filename.str() << "\"" << std::endl;
196 
198  processTreeData<float>(iEvent, file);
199  else
200  processTreeData<double>(iEvent, file);
201 }
202 
203 // ------------ method called once each job just before starting event loop
205 
206 // ------------ method called once each job just after ending the event loop
208 
209 //define this as a plug-in
RunNumber_t run() const
Definition: EventID.h:38
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:40
bool isSparse() const
Definition: PattRecoTree.h:27
Class for storing FFTJet sparse clustering trees.
Definition: PattRecoTree.h:20
fftjet::OpenDXPeakTree< long, fftjet::SparseClusteringTree > SparseFormatter
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
const edm::InputTag treeLabel
#define nullptr
SparseTree sparseTree
std::unique_ptr< fftjet::Functor1< double, fftjet::Peak > > fftjet_PeakFunctor_parser(const edm::ParameterSet &ps)
fftjet::Functor1< double, fftjet::Peak > PeakProperty
unsigned long long EventNumber_t
void checkConfig(const Ptr &ptr, const char *message)
void processTreeData(const edm::Event &, std::ofstream &)
reco::PattRecoTree< Real, reco::PattRecoPeak< Real > > StoredTree
const double etaMax
example_stream void analyze(const edm::Event &, const edm::EventSetup &) override
edm::EDGetTokenT< StoredTree > treeToken
void beginJob()
Definition: Breakpoints.cc:14
static std::string const input
Definition: EdmProvDump.cc:48
fftjet::OpenDXPeakTree< long, fftjet::AbsClusteringTree > DXFormatter
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
ClusteringTree * clusteringTree
~FFTJetTreeDump() override
fftjet::ProximityClusteringTree< fftjet::Peak, long > ClusteringTree
std::unique_ptr< fftjet::AbsDistanceCalculator< fftjet::Peak > > distanceCalc
const bool storeInSinglePrecision
fftjet::SparseClusteringTree< fftjet::Peak, long > SparseTree
const std::string outputPrefix
std::unique_ptr< DXFormatter > denseFormatter
std::unique_ptr< fftjet::AbsDistanceCalculator< fftjet::Peak > > fftjet_DistanceCalculator_parser(const edm::ParameterSet &ps)
std::unique_ptr< PeakProperty > glyphSize
void endJob() override
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)
std::unique_ptr< SparseFormatter > sparseFormatter
edm::EventID id() const
Definition: EventBase.h:59
const double completeEventScale
HLT enums.
const bool insertCompleteEvent
FFTJetTreeDump()=delete
void beginJob() override
unsigned int RunNumber_t
std::unique_ptr< std::vector< double > > iniScales
void analyze(const edm::Event &, const edm::EventSetup &) override
std::unique_ptr< std::vector< double > > fftjet_ScaleSet_parser(const edm::ParameterSet &ps)
std::unique_ptr< PeakProperty > glyphColor