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 // $Id: FFTJetTreeDump.cc,v 1.2 2010/12/07 00:19:43 igv Exp $
17 //
18 //
19 
20 #include <iostream>
21 #include <sstream>
22 #include <fstream>
23 #include <functional>
24 
25 // FFTJet headers
26 #include "fftjet/ProximityClusteringTree.hh"
27 #include "fftjet/OpenDXPeakTree.hh"
28 
29 // user include files
32 
35 
37 
38 // parameter parser header
40 
41 // functions which manipulate storable trees
43 
44 using namespace fftjetcms;
45 
46 //
47 // class declaration
48 //
50 {
51 public:
52  explicit FFTJetTreeDump(const edm::ParameterSet&);
53  ~FFTJetTreeDump();
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;
62 
65  FFTJetTreeDump& operator=(const FFTJetTreeDump&);
66 
67  virtual void beginJob() ;
68  virtual void analyze(const edm::Event&, const edm::EventSetup&);
69  virtual void endJob() ;
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  const double etaMax;
89  const bool insertCompleteEvent;
90  const double completeEventScale;
91 
92  // Distance calculator for the clustering tree
93  std::auto_ptr<fftjet::AbsDistanceCalculator<fftjet::Peak> > distanceCalc;
94 
95  // Scales used
96  std::auto_ptr<std::vector<double> > iniScales;
97 
98  // The sparse clustering tree
100 
101  // Functors which define OpenDX glyph size and color
102  std::auto_ptr<PeakProperty> glyphSize;
103  std::auto_ptr<PeakProperty> glyphColor;
104 
105  // OpenDX formatters
106  std::auto_ptr<DXFormatter> denseFormatter;
107  std::auto_ptr<SparseFormatter> sparseFormatter;
108 
109  unsigned counter;
110 };
111 
112 //
113 // constructors and destructor
114 //
116  : clusteringTree(0),
117  treeLabel(ps.getParameter<edm::InputTag>("treeLabel")),
118  outputPrefix(ps.getParameter<std::string>("outputPrefix")),
119  etaMax(ps.getParameter<double>("etaMax")),
120  storeInSinglePrecision(true),
121  insertCompleteEvent(ps.getParameter<bool>("insertCompleteEvent")),
122  completeEventScale(ps.getParameter<double>("completeEventScale")),
123  counter(0)
124 {
125  if (etaMax < 0.0)
126  throw cms::Exception("FFTJetBadConfig")
127  << "etaMax can not be negative" << std::endl;
128 
129  // Build the set of pattern recognition scales
130  const edm::ParameterSet& InitialScales(
131  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(
138  ps.getParameter<edm::ParameterSet>("TreeDistanceCalculator"));
139  distanceCalc = fftjet_DistanceCalculator_parser(TreeDistanceCalculator);
140  checkConfig(distanceCalc, "invalid tree distance calculator");
141 
142  // Determine representations for the OpenDX glyph size and color
143  const edm::ParameterSet& GlyphSize(
144  ps.getParameter<edm::ParameterSet>("GlyphSize"));
146  checkConfig(glyphSize, "invalid glyph size parameters");
147 
148  const edm::ParameterSet& GlyphColor(
149  ps.getParameter<edm::ParameterSet>("GlyphColor"));
151  checkConfig(glyphColor, "invalid glyph color parameters");
152 
153  // Build the tree formatters
154  denseFormatter = std::auto_ptr<DXFormatter>(new DXFormatter(
155  glyphSize.get(), glyphColor.get(), etaMax));
156  sparseFormatter = std::auto_ptr<SparseFormatter>(new SparseFormatter(
157  glyphSize.get(), glyphColor.get(), etaMax));
158 
159  // Build the clustering tree
161 }
162 
163 
165 {
166  delete clusteringTree;
167 }
168 
169 
170 //
171 // member functions
172 //
173 template<class Real>
175  std::ofstream& file)
176 {
178 
179  // Get the event number
180  const unsigned long runNum = iEvent.id().run();
181  const unsigned long evNum = iEvent.id().event();
182 
183  // Get the input
185  iEvent.getByLabel(treeLabel, input);
186 
187  const double eventScale = insertCompleteEvent ? completeEventScale : 0.0;
188  if (input->isSparse())
189  {
191  eventScale, &sparseTree);
192  sparseFormatter->setTree(sparseTree, runNum, evNum);
193  file << *sparseFormatter << std::endl;
194  }
195  else
196  {
197  densePeakTreeFromStorable(*input, iniScales.get(),
198  eventScale, clusteringTree);
199  denseFormatter->setTree(*clusteringTree, runNum, evNum);
200  file << *denseFormatter << std::endl;
201  }
202 }
203 
204 
205 // ------------ method called to for each event ------------
207  const edm::EventSetup& iSetup)
208 {
209  // Create the file name
210  std::ostringstream filename;
211  filename << outputPrefix << '_' << counter++ << ".dx";
212 
213  // Open the file
214  std::ofstream file(filename.str().c_str());
215  if (!file)
216  throw cms::Exception("FFTJetBadConfig")
217  << "Failed to open file \""
218  << filename.str() << "\"" << std::endl;
219 
221  processTreeData<float>(iEvent, file);
222  else
223  processTreeData<double>(iEvent, file);
224 }
225 
226 
227 // ------------ method called once each job just before starting event loop
229 {
230 }
231 
232 // ------------ method called once each job just after ending the event loop
234 }
235 
236 //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:19
virtual void endJob()
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
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
virtual void beginJob()
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:361
fftjet::OpenDXPeakTree< long, fftjet::SparseClusteringTree > SparseFormatter
std::auto_ptr< fftjet::Functor1< double, fftjet::Peak > > fftjet_PeakFunctor_parser(const edm::ParameterSet &ps)
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)
const bool insertCompleteEvent
tuple filename
Definition: lut2db_cfg.py:20
virtual void analyze(const edm::Event &, const edm::EventSetup &)
std::auto_ptr< std::vector< double > > fftjet_ScaleSet_parser(const edm::ParameterSet &ps)