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 {
52 public:
53  explicit FFTJetTreeDump(const edm::ParameterSet&);
54  ~FFTJetTreeDump() override;
55 
56 private:
57  // Useful local typedefs
58  typedef fftjet::ProximityClusteringTree<fftjet::Peak,long> ClusteringTree;
59  typedef fftjet::SparseClusteringTree<fftjet::Peak,long> SparseTree;
60  typedef fftjet::OpenDXPeakTree<long,fftjet::AbsClusteringTree> DXFormatter;
61  typedef fftjet::OpenDXPeakTree<long,fftjet::SparseClusteringTree> SparseFormatter;
62  typedef fftjet::Functor1<double,fftjet::Peak> PeakProperty;
64 
65  FFTJetTreeDump() = delete;
66  FFTJetTreeDump(const FFTJetTreeDump&) = delete;
67  FFTJetTreeDump& operator=(const FFTJetTreeDump&) = delete;
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  {
79  if (ptr.get() == nullptr)
80  throw cms::Exception("FFTJetBadConfig") << message << std::endl;
81  }
82 
83  // ----------member data ---------------------------
84  // The complete clustering tree
85  ClusteringTree* clusteringTree;
86 
89 
91  const double etaMax;
93  const bool insertCompleteEvent;
94  const double completeEventScale;
95 
96  // Distance calculator for the clustering tree
97  std::auto_ptr<fftjet::AbsDistanceCalculator<fftjet::Peak> > distanceCalc;
98 
99  // Scales used
100  std::auto_ptr<std::vector<double> > iniScales;
101 
102  // The sparse clustering tree
103  SparseTree sparseTree;
104 
105  // Functors which define OpenDX glyph size and color
106  std::auto_ptr<PeakProperty> glyphSize;
107  std::auto_ptr<PeakProperty> glyphColor;
108 
109  // OpenDX formatters
110  std::auto_ptr<DXFormatter> denseFormatter;
111  std::auto_ptr<SparseFormatter> sparseFormatter;
112 
113  unsigned counter;
114 };
115 
116 //
117 // constructors and destructor
118 //
120  : clusteringTree(nullptr),
121  treeLabel(ps.getParameter<edm::InputTag>("treeLabel")),
122  outputPrefix(ps.getParameter<std::string>("outputPrefix")),
123  etaMax(ps.getParameter<double>("etaMax")),
124  storeInSinglePrecision(true),
125  insertCompleteEvent(ps.getParameter<bool>("insertCompleteEvent")),
126  completeEventScale(ps.getParameter<double>("completeEventScale")),
127  counter(0)
128 {
129  if (etaMax < 0.0)
130  throw cms::Exception("FFTJetBadConfig")
131  << "etaMax can not be negative" << std::endl;
132 
133  // Build the set of pattern recognition scales
135  ps.getParameter<edm::ParameterSet>("InitialScales"));
136  iniScales = fftjet_ScaleSet_parser(InitialScales);
137  checkConfig(iniScales, "invalid set of scales");
138  std::sort(iniScales->begin(), iniScales->end(), std::greater<double>());
139 
140  // Distance calculator for the clustering tree
142  ps.getParameter<edm::ParameterSet>("TreeDistanceCalculator"));
143  distanceCalc = fftjet_DistanceCalculator_parser(TreeDistanceCalculator);
144  checkConfig(distanceCalc, "invalid tree distance calculator");
145 
146  // Determine representations for the OpenDX glyph size and color
148  ps.getParameter<edm::ParameterSet>("GlyphSize"));
150  checkConfig(glyphSize, "invalid glyph size parameters");
151 
153  ps.getParameter<edm::ParameterSet>("GlyphColor"));
155  checkConfig(glyphColor, "invalid glyph color parameters");
156 
157  // Build the tree formatters
158  denseFormatter = std::auto_ptr<DXFormatter>(new DXFormatter(
159  glyphSize.get(), glyphColor.get(), etaMax));
160  sparseFormatter = std::auto_ptr<SparseFormatter>(new SparseFormatter(
161  glyphSize.get(), glyphColor.get(), etaMax));
162 
163  // Build the clustering tree
165 
166  treeToken = consumes<StoredTree>(treeLabel);
167 }
168 
169 
171 {
172  delete clusteringTree;
173 }
174 
175 
176 //
177 // member functions
178 //
179 template<class Real>
181  std::ofstream& file)
182 {
183  // Get the event number
184  edm::RunNumber_t const runNum = iEvent.id().run();
185  edm::EventNumber_t const evNum = iEvent.id().event();
186 
187  // Get the input
189  iEvent.getByToken(treeToken, input);
190 
191  const double eventScale = insertCompleteEvent ? completeEventScale : 0.0;
192  if (input->isSparse())
193  {
195  eventScale, &sparseTree);
196  sparseFormatter->setTree(sparseTree, runNum, evNum);
197  file << *sparseFormatter << std::endl;
198  }
199  else
200  {
201  densePeakTreeFromStorable(*input, iniScales.get(),
202  eventScale, clusteringTree);
203  denseFormatter->setTree(*clusteringTree, runNum, evNum);
204  file << *denseFormatter << std::endl;
205  }
206 }
207 
208 
209 // ------------ method called to for each event ------------
211  const edm::EventSetup& iSetup)
212 {
213  // Create the file name
214  std::ostringstream filename;
215  filename << outputPrefix << '_' << counter++ << ".dx";
216 
217  // Open the file
218  std::ofstream file(filename.str().c_str());
219  if (!file)
220  throw cms::Exception("FFTJetBadConfig")
221  << "Failed to open file \""
222  << filename.str() << "\"" << std::endl;
223 
225  processTreeData<float>(iEvent, file);
226  else
227  processTreeData<double>(iEvent, file);
228 }
229 
230 
231 // ------------ method called once each job just before starting event loop
233 {
234 }
235 
236 // ------------ method called once each job just after ending the event loop
238 }
239 
240 //define this as a plug-in
RunNumber_t run() const
Definition: EventID.h:39
fftjet::OpenDXPeakTree< long, fftjet::AbsClusteringTree > DXFormatter
std::auto_ptr< PeakProperty > glyphSize
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:41
bool isSparse() const
Definition: PattRecoTree.h:28
Class for storing FFTJet sparse clustering trees.
Definition: PattRecoTree.h:20
std::auto_ptr< SparseFormatter > sparseFormatter
virtual example_stream void analyze(const edm::Event &, const edm::EventSetup &) override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
const edm::InputTag treeLabel
SparseTree sparseTree
unsigned long long EventNumber_t
void checkConfig(const Ptr &ptr, const char *message)
void processTreeData(const edm::Event &, std::ofstream &)
#define nullptr
const double etaMax
edm::EDGetTokenT< StoredTree > treeToken
void beginJob()
Definition: Breakpoints.cc:15
static std::string const input
Definition: EdmProvDump.cc:44
std::auto_ptr< std::vector< double > > iniScales
std::auto_ptr< DXFormatter > denseFormatter
int iEvent
Definition: GenABIO.cc:230
ClusteringTree * clusteringTree
fftjet::Functor1< double, fftjet::Peak > PeakProperty
~FFTJetTreeDump() override
fftjet::ProximityClusteringTree< fftjet::Peak, long > ClusteringTree
std::auto_ptr< fftjet::AbsDistanceCalculator< fftjet::Peak > > distanceCalc
reco::PattRecoTree< Real, reco::PattRecoPeak< Real > > StoredTree
const bool storeInSinglePrecision
std::auto_ptr< PeakProperty > glyphColor
const std::string outputPrefix
fftjet::SparseClusteringTree< fftjet::Peak, long > SparseTree
fftjet::OpenDXPeakTree< long, fftjet::SparseClusteringTree > SparseFormatter
std::auto_ptr< fftjet::Functor1< double, fftjet::Peak > > fftjet_PeakFunctor_parser(const edm::ParameterSet &ps)
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)
edm::EventID id() const
Definition: EventBase.h:60
const double completeEventScale
HLT enums.
std::auto_ptr< fftjet::AbsDistanceCalculator< fftjet::Peak > > fftjet_DistanceCalculator_parser(const edm::ParameterSet &ps)
const bool insertCompleteEvent
FFTJetTreeDump()=delete
void beginJob() override
unsigned int RunNumber_t
void analyze(const edm::Event &, const edm::EventSetup &) override
std::auto_ptr< std::vector< double > > fftjet_ScaleSet_parser(const edm::ParameterSet &ps)