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;
62 
65  FFTJetTreeDump& operator=(const FFTJetTreeDump&);
66 
67  virtual void beginJob() override ;
68  virtual void analyze(const edm::Event&, const edm::EventSetup&) override;
69  virtual void endJob() override ;
70 
71  template<class Real>
72  void processTreeData(const edm::Event&, std::ofstream&);
73 
74  template<class Ptr>
75  void checkConfig(const Ptr& ptr, const char* message)
76  {
77  if (ptr.get() == NULL)
78  throw cms::Exception("FFTJetBadConfig") << message << std::endl;
79  }
80 
81  // ----------member data ---------------------------
82  // The complete clustering tree
84 
87 
89  const double etaMax;
91  const bool insertCompleteEvent;
92  const double completeEventScale;
93 
94  // Distance calculator for the clustering tree
95  std::auto_ptr<fftjet::AbsDistanceCalculator<fftjet::Peak> > distanceCalc;
96 
97  // Scales used
98  std::auto_ptr<std::vector<double> > iniScales;
99 
100  // The sparse clustering tree
102 
103  // Functors which define OpenDX glyph size and color
104  std::auto_ptr<PeakProperty> glyphSize;
105  std::auto_ptr<PeakProperty> glyphColor;
106 
107  // OpenDX formatters
108  std::auto_ptr<DXFormatter> denseFormatter;
109  std::auto_ptr<SparseFormatter> sparseFormatter;
110 
111  unsigned counter;
112 };
113 
114 //
115 // constructors and destructor
116 //
118  : clusteringTree(0),
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 {
127  if (etaMax < 0.0)
128  throw cms::Exception("FFTJetBadConfig")
129  << "etaMax can not be negative" << std::endl;
130 
131  // Build the set of pattern recognition scales
132  const edm::ParameterSet& InitialScales(
133  ps.getParameter<edm::ParameterSet>("InitialScales"));
134  iniScales = fftjet_ScaleSet_parser(InitialScales);
135  checkConfig(iniScales, "invalid set of scales");
136  std::sort(iniScales->begin(), iniScales->end(), std::greater<double>());
137 
138  // Distance calculator for the clustering tree
139  const edm::ParameterSet& TreeDistanceCalculator(
140  ps.getParameter<edm::ParameterSet>("TreeDistanceCalculator"));
141  distanceCalc = fftjet_DistanceCalculator_parser(TreeDistanceCalculator);
142  checkConfig(distanceCalc, "invalid tree distance calculator");
143 
144  // Determine representations for the OpenDX glyph size and color
145  const edm::ParameterSet& GlyphSize(
146  ps.getParameter<edm::ParameterSet>("GlyphSize"));
148  checkConfig(glyphSize, "invalid glyph size parameters");
149 
150  const edm::ParameterSet& GlyphColor(
151  ps.getParameter<edm::ParameterSet>("GlyphColor"));
153  checkConfig(glyphColor, "invalid glyph color parameters");
154 
155  // Build the tree formatters
156  denseFormatter = std::auto_ptr<DXFormatter>(new DXFormatter(
157  glyphSize.get(), glyphColor.get(), etaMax));
158  sparseFormatter = std::auto_ptr<SparseFormatter>(new SparseFormatter(
159  glyphSize.get(), glyphColor.get(), etaMax));
160 
161  // Build the clustering tree
163 
164  treeToken = consumes<StoredTree>(treeLabel);
165 }
166 
167 
169 {
170  delete clusteringTree;
171 }
172 
173 
174 //
175 // member functions
176 //
177 template<class Real>
179  std::ofstream& file)
180 {
181  // Get the event number
182  const unsigned long runNum = iEvent.id().run();
183  const unsigned long evNum = iEvent.id().event();
184 
185  // Get the input
187  iEvent.getByToken(treeToken, input);
188 
189  const double eventScale = insertCompleteEvent ? completeEventScale : 0.0;
190  if (input->isSparse())
191  {
193  eventScale, &sparseTree);
194  sparseFormatter->setTree(sparseTree, runNum, evNum);
195  file << *sparseFormatter << std::endl;
196  }
197  else
198  {
199  densePeakTreeFromStorable(*input, iniScales.get(),
200  eventScale, clusteringTree);
201  denseFormatter->setTree(*clusteringTree, runNum, evNum);
202  file << *denseFormatter << std::endl;
203  }
204 }
205 
206 
207 // ------------ method called to for each event ------------
209  const edm::EventSetup& iSetup)
210 {
211  // Create the file name
212  std::ostringstream filename;
213  filename << outputPrefix << '_' << counter++ << ".dx";
214 
215  // Open the file
216  std::ofstream file(filename.str().c_str());
217  if (!file)
218  throw cms::Exception("FFTJetBadConfig")
219  << "Failed to open file \""
220  << filename.str() << "\"" << std::endl;
221 
223  processTreeData<float>(iEvent, file);
224  else
225  processTreeData<double>(iEvent, file);
226 }
227 
228 
229 // ------------ method called once each job just before starting event loop
231 {
232 }
233 
234 // ------------ method called once each job just after ending the event loop
236 }
237 
238 //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
Class for storing FFTJet sparse clustering trees.
Definition: PattRecoTree.h:18
std::auto_ptr< SparseFormatter > sparseFormatter
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
#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
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
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)
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)