CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FFTJetCorrectionESProducer.h
Go to the documentation of this file.
1 #ifndef JetMETCorrections_FFTJetModules_plugins_FFTJetCorrectionESProducer_h
2 #define JetMETCorrections_FFTJetModules_plugins_FFTJetCorrectionESProducer_h
3 
4 // -*- C++ -*-
5 //
6 // Package: FFTJetModules
7 // Class: FFTJetCorrectionESProducer
8 //
16 //
17 // Original Author: Igor Volobouev
18 // Created: Thu Aug 2 22:34:02 CDT 2012
19 //
20 //
21 
22 
23 // system include files
24 #include <sstream>
25 
26 #include "boost/shared_ptr.hpp"
27 
28 #include "Alignment/Geners/interface/CompressedIO.hh"
29 
30 // user include files
34 
37 
42 
44 
45 //
46 // Don't make the function which builds the sequence
47 // a member of the producer class -- this way it gets
48 // instantiated only as many times as there are corrector
49 // sequence types but not as many times as there are
50 // record types.
51 //
52 template<class CorrectorSequence>
53 static boost::shared_ptr<CorrectorSequence>
55  const FFTJetCorrectorParameters& tablePars,
56  const std::vector<edm::ParameterSet>& sequence,
57  const bool isArchiveCompressed, const bool verbose)
58 {
59  typedef typename CorrectorSequence::Corrector Corrector;
60  typedef typename CorrectorSequence::jet_type jet_type;
61 
62  // Load the archive stored in the FFTJetCorrectorParameters object
63  CPP11_auto_ptr<gs::StringArchive> ar;
64  {
65  std::istringstream is(tablePars.str());
66  if (isArchiveCompressed)
67  ar = gs::read_compressed_item<gs::StringArchive>(is);
68  else
69  ar = gs::read_item<gs::StringArchive>(is);
70  }
71 
72  // Create an empty corrector sequence
73  boost::shared_ptr<CorrectorSequence> ptr(new CorrectorSequence());
74 
75  // Go over the parameter sets in the VPSet and
76  // configure all correction levels. Add each new
77  // level to the sequence.
78  const unsigned nLevels = sequence.size();
79  for (unsigned lev=0; lev<nLevels; ++lev)
80  ptr->addCorrector(parseFFTJetCorrector<Corrector>(
81  sequence[lev], *ar, verbose));
82 
83  // Check for proper level order.
84  // Assume that level 0 is special and can happen anywhere.
85  unsigned previousLevel = 0;
86  for (unsigned lev=0; lev<nLevels; ++lev)
87  {
88  const unsigned level = (*ptr)[lev].level();
89  if (level)
90  {
91  if (level <= previousLevel)
92  throw cms::Exception("FFTJetBadConfig")
93  << "Error in buildCorrectorSequence: "
94  << "correction levels are out of order\n";
95  previousLevel = level;
96  }
97  }
98 
99  return ptr;
100 }
101 
102 //
103 // class declaration
104 //
105 template<typename CT>
107 {
108 public:
109  // Lookup various types
111  typedef FFTJetCorrectorSequence<
112  jet_type,
116  typedef boost::shared_ptr<CorrectorSequence> ReturnType;
119 
122 
123  ReturnType produce(const MyRecord&);
124 
125 private:
126  inline void doWhenChanged(const ParentRecord&)
127  {remakeProduct = true;}
128 
129  // Module parameters
130  std::vector<edm::ParameterSet> sequence;
132  bool verbose;
133 
134  // Other module variables
137 };
138 
139 //
140 // constructors and destructor
141 //
142 template<typename CT>
144  const edm::ParameterSet& psIn)
145  : sequence(psIn.getParameter<std::vector<edm::ParameterSet> >("sequence")),
146  isArchiveCompressed(psIn.getParameter<bool>("isArchiveCompressed")),
147  verbose(psIn.getUntrackedParameter<bool>("verbose", false)),
148  remakeProduct(true)
149 {
150  // The following line is needed to tell the framework what
151  // data is being produced
153 }
154 
155 // ------------ method called to produce the data ------------
156 template<typename CT>
159 {
160  if (remakeProduct)
161  {
162  // According to:
163  // https://twiki.cern.ch/twiki/bin/view/CMSPublic/SWGuideHowToGetDependentRecord
164  //
165  // If ARecord is dependent on BRecord then you can
166  // call the method getRecord of ARecord:
167  //
168  // const BRecord& b = aRecord.getRecord<BRecord>();
169  //
170  const ParentRecord& rec = iRecord.template getRecord<ParentRecord>();
172  rec.get(parHandle);
173  product = buildCorrectorSequence<CorrectorSequence>(
174  *parHandle, sequence, isArchiveCompressed, verbose);
175  remakeProduct = false;
176  }
177  return product;
178 }
179 
180 #endif // JetMETCorrections_FFTJetModules_plugins_FFTJetCorrectionESProducer_h
boost::shared_ptr< CorrectorSequence > ReturnType
depends_on::OneHolder< T, TDependsOnRecord > dependsOn(void(T::*iT)(const TDependsOnRecord &))
void doWhenChanged(const ParentRecord &)
bool verbose
const std::string & str() const
void setWhatProduced(T *iThis, const es::Label &iLabel=es::Label())
Definition: ESProducer.h:115
FFTJetCorrectorParametersRcd< CT > ParentRecord
FFTJetCorrectionESProducer(const edm::ParameterSet &)
void get(HolderT &iHolder) const
FFTJetCorrectionsTypemap< CT >::jet_type jet_type
ReturnType produce(const MyRecord &)
FFTJetCorrectorSequenceRcd< CT > MyRecord
static boost::shared_ptr< CorrectorSequence > buildCorrectorSequence(const FFTJetCorrectorParameters &tablePars, const std::vector< edm::ParameterSet > &sequence, const bool isArchiveCompressed, const bool verbose)
tuple level
Definition: testEve_cfg.py:34
volatile std::atomic< bool > shutdown_flag false
FFTJetCorrectorSequence< jet_type, FFTJetCorrectorTransientFromJet, FFTJetCorrectorResultFromTransient > CorrectorSequence
std::vector< edm::ParameterSet > sequence