CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros 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 <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 
43 using namespace fftjetcms;
44 
45 //
46 // class declaration
47 //
49 {
50 public:
51  explicit FFTJetTreeDump(const edm::ParameterSet&);
52  ~FFTJetTreeDump();
53 
54 private:
55  // Useful local typedefs
56  typedef fftjet::ProximityClusteringTree<fftjet::Peak,long> ClusteringTree;
57  typedef fftjet::SparseClusteringTree<fftjet::Peak,long> SparseTree;
58  typedef fftjet::OpenDXPeakTree<long,fftjet::AbsClusteringTree> DXFormatter;
59  typedef fftjet::OpenDXPeakTree<long,fftjet::SparseClusteringTree> SparseFormatter;
60  typedef fftjet::Functor1<double,fftjet::Peak> PeakProperty;
61 
64  FFTJetTreeDump& operator=(const FFTJetTreeDump&);
65 
66  virtual void beginJob() override ;
67  virtual void analyze(const edm::Event&, const edm::EventSetup&) override;
68  virtual void endJob() override ;
69 
70  template<class Real>
71  void processTreeData(const edm::Event&, std::ofstream&);
72 
73  template<class Ptr>
74  void checkConfig(const Ptr& ptr, const char* message)
75  {
76  if (ptr.get() == NULL)
77  throw cms::Exception("FFTJetBadConfig") << message << std::endl;
78  }
79 
80  // ----------member data ---------------------------
81  // The complete clustering tree
83 
86  const double etaMax;
88  const bool insertCompleteEvent;
89  const double completeEventScale;
90 
91  // Distance calculator for the clustering tree
92  std::auto_ptr<fftjet::AbsDistanceCalculator<fftjet::Peak> > distanceCalc;
93 
94  // Scales used
95  std::auto_ptr<std::vector<double> > iniScales;
96 
97  // The sparse clustering tree
99 
100  // Functors which define OpenDX glyph size and color
101  std::auto_ptr<PeakProperty> glyphSize;
102  std::auto_ptr<PeakProperty> glyphColor;
103 
104  // OpenDX formatters
105  std::auto_ptr<DXFormatter> denseFormatter;
106  std::auto_ptr<SparseFormatter> sparseFormatter;
107 
108  unsigned counter;
109 };
110 
111 //
112 // constructors and destructor
113 //
115  : clusteringTree(0),
116  treeLabel(ps.getParameter<edm::InputTag>("treeLabel")),
117  outputPrefix(ps.getParameter<std::string>("outputPrefix")),
118  etaMax(ps.getParameter<double>("etaMax")),
119  storeInSinglePrecision(true),
120  insertCompleteEvent(ps.getParameter<bool>("insertCompleteEvent")),
121  completeEventScale(ps.getParameter<double>("completeEventScale")),
122  counter(0)
123 {
124  if (etaMax < 0.0)
125  throw cms::Exception("FFTJetBadConfig")
126  << "etaMax can not be negative" << std::endl;
127 
128  // Build the set of pattern recognition scales
129  const edm::ParameterSet& InitialScales(
130  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(
137  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(
143  ps.getParameter<edm::ParameterSet>("GlyphSize"));
145  checkConfig(glyphSize, "invalid glyph size parameters");
146 
147  const edm::ParameterSet& GlyphColor(
148  ps.getParameter<edm::ParameterSet>("GlyphColor"));
150  checkConfig(glyphColor, "invalid glyph color parameters");
151 
152  // Build the tree formatters
153  denseFormatter = std::auto_ptr<DXFormatter>(new DXFormatter(
154  glyphSize.get(), glyphColor.get(), etaMax));
155  sparseFormatter = std::auto_ptr<SparseFormatter>(new SparseFormatter(
156  glyphSize.get(), glyphColor.get(), etaMax));
157 
158  // Build the clustering tree
160 }
161 
162 
164 {
165  delete clusteringTree;
166 }
167 
168 
169 //
170 // member functions
171 //
172 template<class Real>
174  std::ofstream& file)
175 {
177 
178  // Get the event number
179  const unsigned long runNum = iEvent.id().run();
180  const unsigned long evNum = iEvent.id().event();
181 
182  // Get the input
184  iEvent.getByLabel(treeLabel, input);
185 
186  const double eventScale = insertCompleteEvent ? completeEventScale : 0.0;
187  if (input->isSparse())
188  {
190  eventScale, &sparseTree);
191  sparseFormatter->setTree(sparseTree, runNum, evNum);
192  file << *sparseFormatter << std::endl;
193  }
194  else
195  {
196  densePeakTreeFromStorable(*input, iniScales.get(),
197  eventScale, clusteringTree);
198  denseFormatter->setTree(*clusteringTree, runNum, evNum);
199  file << *denseFormatter << std::endl;
200  }
201 }
202 
203 
204 // ------------ method called to for each event ------------
206  const edm::EventSetup& iSetup)
207 {
208  // Create the file name
209  std::ostringstream filename;
210  filename << outputPrefix << '_' << counter++ << ".dx";
211 
212  // Open the file
213  std::ofstream file(filename.str().c_str());
214  if (!file)
215  throw cms::Exception("FFTJetBadConfig")
216  << "Failed to open file \""
217  << filename.str() << "\"" << std::endl;
218 
220  processTreeData<float>(iEvent, file);
221  else
222  processTreeData<double>(iEvent, file);
223 }
224 
225 
226 // ------------ method called once each job just before starting event loop
228 {
229 }
230 
231 // ------------ method called once each job just after ending the event loop
233 }
234 
235 //define this as a plug-in
RunNumber_t run() const
Definition: EventID.h:42
fftjet::OpenDXPeakTree< long, fftjet::AbsClusteringTree > DXFormatter
std::auto_ptr< PeakProperty > glyphSize
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:44
Class for storing FFTJet sparse clustering trees.
Definition: PattRecoTree.h:18
std::auto_ptr< SparseFormatter > sparseFormatter
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
const edm::InputTag treeLabel
SparseTree sparseTree
void checkConfig(const Ptr &ptr, const char *message)
#define NULL
Definition: scimark2.h:8
void processTreeData(const edm::Event &, std::ofstream &)
const double etaMax
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:243
ClusteringTree * clusteringTree
fftjet::Functor1< double, fftjet::Peak > PeakProperty
fftjet::ProximityClusteringTree< fftjet::Peak, long > ClusteringTree
std::auto_ptr< fftjet::AbsDistanceCalculator< fftjet::Peak > > distanceCalc
const bool storeInSinglePrecision
std::auto_ptr< PeakProperty > glyphColor
const std::string outputPrefix
fftjet::SparseClusteringTree< fftjet::Peak, long > SparseTree
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
fftjet::OpenDXPeakTree< long, fftjet::SparseClusteringTree > SparseFormatter
std::auto_ptr< fftjet::Functor1< double, fftjet::Peak > > fftjet_PeakFunctor_parser(const edm::ParameterSet &ps)
virtual 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:56
const double completeEventScale
std::auto_ptr< fftjet::AbsDistanceCalculator< fftjet::Peak > > fftjet_DistanceCalculator_parser(const edm::ParameterSet &ps)
static std::atomic< unsigned int > counter
const bool insertCompleteEvent
tuple filename
Definition: lut2db_cfg.py:20
virtual void beginJob() override
virtual void analyze(const edm::Event &, const edm::EventSetup &) override
std::auto_ptr< std::vector< double > > fftjet_ScaleSet_parser(const edm::ParameterSet &ps)