CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 <memory>
21 
22 #include <fstream>
23 #include <functional>
24 #include <sstream>
25 
26 // FFTJet headers
27 #include "fftjet/ProximityClusteringTree.hh"
28 #include "fftjet/OpenDXPeakTree.hh"
29 
30 // user include files
33 
36 
38 
39 // parameter parser header
41 
42 // functions which manipulate storable trees
44 
46 
47 using namespace fftjetcms;
48 
49 //
50 // class declaration
51 //
53 public:
54  explicit FFTJetTreeDump(const edm::ParameterSet&);
55  FFTJetTreeDump() = delete;
56  FFTJetTreeDump(const FFTJetTreeDump&) = delete;
57  FFTJetTreeDump& operator=(const FFTJetTreeDump&) = delete;
58  ~FFTJetTreeDump() override;
59 
60 private:
61  // Useful local typedefs
62  typedef fftjet::ProximityClusteringTree<fftjet::Peak, long> ClusteringTree;
63  typedef fftjet::SparseClusteringTree<fftjet::Peak, long> SparseTree;
64  typedef fftjet::OpenDXPeakTree<long, fftjet::AbsClusteringTree> DXFormatter;
65  typedef fftjet::OpenDXPeakTree<long, fftjet::SparseClusteringTree> SparseFormatter;
66  typedef fftjet::Functor1<double, fftjet::Peak> PeakProperty;
68 
69  void beginJob() override;
70  void analyze(const edm::Event&, const edm::EventSetup&) override;
71  void endJob() override;
72 
73  template <class Real>
74  void processTreeData(const edm::Event&, std::ofstream&);
75 
76  template <class Ptr>
77  void checkConfig(const Ptr& ptr, const char* message) {
78  if (ptr.get() == nullptr)
79  throw cms::Exception("FFTJetBadConfig") << message << std::endl;
80  }
81 
82  // ----------member data ---------------------------
83  // The complete clustering tree
85 
88 
90  const double etaMax;
92  const bool insertCompleteEvent;
93  const double completeEventScale;
94 
95  // Distance calculator for the clustering tree
96  std::unique_ptr<fftjet::AbsDistanceCalculator<fftjet::Peak> > distanceCalc;
97 
98  // Scales used
99  std::unique_ptr<std::vector<double> > iniScales;
100 
101  // The sparse clustering tree
103 
104  // Functors which define OpenDX glyph size and color
105  std::unique_ptr<PeakProperty> glyphSize;
106  std::unique_ptr<PeakProperty> glyphColor;
107 
108  // OpenDX formatters
109  std::unique_ptr<DXFormatter> denseFormatter;
110  std::unique_ptr<SparseFormatter> sparseFormatter;
111 
112  unsigned counter;
113 };
114 
115 //
116 // constructors and destructor
117 //
119  : clusteringTree(nullptr),
120  treeLabel(ps.getParameter<edm::InputTag>("treeLabel")),
121  outputPrefix(ps.getParameter<std::string>("outputPrefix")),
122  etaMax(ps.getParameter<double>("etaMax")),
123  storeInSinglePrecision(true),
124  insertCompleteEvent(ps.getParameter<bool>("insertCompleteEvent")),
125  completeEventScale(ps.getParameter<double>("completeEventScale")),
126  counter(0) {
127  if (etaMax < 0.0)
128  throw cms::Exception("FFTJetBadConfig") << "etaMax can not be negative" << std::endl;
129 
130  // Build the set of pattern recognition scales
131  const edm::ParameterSet& InitialScales(ps.getParameter<edm::ParameterSet>("InitialScales"));
132  iniScales = fftjet_ScaleSet_parser(InitialScales);
133  checkConfig(iniScales, "invalid set of scales");
134  std::sort(iniScales->begin(), iniScales->end(), std::greater<double>());
135 
136  // Distance calculator for the clustering tree
137  const edm::ParameterSet& TreeDistanceCalculator(ps.getParameter<edm::ParameterSet>("TreeDistanceCalculator"));
138  distanceCalc = fftjet_DistanceCalculator_parser(TreeDistanceCalculator);
139  checkConfig(distanceCalc, "invalid tree distance calculator");
140 
141  // Determine representations for the OpenDX glyph size and color
142  const edm::ParameterSet& GlyphSize(ps.getParameter<edm::ParameterSet>("GlyphSize"));
144  checkConfig(glyphSize, "invalid glyph size parameters");
145 
146  const edm::ParameterSet& GlyphColor(ps.getParameter<edm::ParameterSet>("GlyphColor"));
148  checkConfig(glyphColor, "invalid glyph color parameters");
149 
150  // Build the tree formatters
151  denseFormatter = std::make_unique<DXFormatter>(glyphSize.get(), glyphColor.get(), etaMax);
152  sparseFormatter = std::make_unique<SparseFormatter>(glyphSize.get(), glyphColor.get(), etaMax);
153 
154  // Build the clustering tree
156 
157  treeToken = consumes<StoredTree>(treeLabel);
158 }
159 
161 
162 //
163 // member functions
164 //
165 template <class Real>
167  // Get the event number
168  edm::RunNumber_t const runNum = iEvent.id().run();
169  edm::EventNumber_t const evNum = iEvent.id().event();
170 
171  // Get the input
173  iEvent.getByToken(treeToken, input);
174 
175  const double eventScale = insertCompleteEvent ? completeEventScale : 0.0;
176  if (input->isSparse()) {
177  sparsePeakTreeFromStorable(*input, iniScales.get(), eventScale, &sparseTree);
178  sparseFormatter->setTree(sparseTree, runNum, evNum);
179  file << *sparseFormatter << std::endl;
180  } else {
181  densePeakTreeFromStorable(*input, iniScales.get(), eventScale, clusteringTree);
182  denseFormatter->setTree(*clusteringTree, runNum, evNum);
183  file << *denseFormatter << std::endl;
184  }
185 }
186 
187 // ------------ method called to for each event ------------
189  // Create the file name
190  std::ostringstream filename;
191  filename << outputPrefix << '_' << counter++ << ".dx";
192 
193  // Open the file
194  std::ofstream file(filename.str().c_str());
195  if (!file)
196  throw cms::Exception("FFTJetBadConfig") << "Failed to open file \"" << filename.str() << "\"" << std::endl;
197 
199  processTreeData<float>(iEvent, file);
200  else
201  processTreeData<double>(iEvent, file);
202 }
203 
204 // ------------ method called once each job just before starting event loop
206 
207 // ------------ method called once each job just after ending the event loop
209 
210 //define this as a plug-in
RunNumber_t run() const
Definition: EventID.h:38
EventNumber_t event() const
Definition: EventID.h:40
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:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const edm::InputTag treeLabel
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:47
fftjet::OpenDXPeakTree< long, fftjet::AbsClusteringTree > DXFormatter
int iEvent
Definition: GenABIO.cc:224
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)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::unique_ptr< SparseFormatter > sparseFormatter
edm::EventID id() const
Definition: EventBase.h:59
const double completeEventScale
static std::atomic< unsigned int > counter
const bool insertCompleteEvent
tuple filename
Definition: lut2db_cfg.py:20
FFTJetTreeDump()=delete
void beginJob() override
unsigned int RunNumber_t
tuple etaMax
Definition: Puppi_cff.py:46
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