CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
JetCorrectionProducer.h
Go to the documentation of this file.
1 #ifndef JetCorrectionProducer_h
2 #define JetCorrectionProducer_h
3 
4 #include <sstream>
5 #include <string>
6 #include <vector>
7 
20 
21 namespace edm
22 {
23  class ParameterSet;
24 }
25 
26 class JetCorrector;
27 
28 namespace cms
29 {
30  template<class T>
32  {
33  public:
34  typedef std::vector<T> JetCollection;
35  explicit JetCorrectionProducer (const edm::ParameterSet& fParameters);
36  virtual ~JetCorrectionProducer () {}
37  virtual void produce(edm::Event&, const edm::EventSetup&);
38  private:
40  std::vector <std::string> mCorrectorNames;
41  // cache
42  std::vector <const JetCorrector*> mCorrectors;
43  unsigned long long mCacheId;
44  bool mVerbose;
45  };
46 }
47 
48 // ---------- implementation ----------
49 
50 namespace cms {
51 
52  template<class T>
54  : mInput(consumes<JetCollection>(fConfig.getParameter <edm::InputTag> ("src")))
55  , mCorrectorNames(fConfig.getParameter<std::vector<std::string> >("correctors"))
56  , mCorrectors(mCorrectorNames.size(), 0)
57  , mCacheId (0)
58  , mVerbose (fConfig.getUntrackedParameter <bool> ("verbose", false))
59  {
60  std::string alias = fConfig.getUntrackedParameter <std::string> ("alias", "");
61  if (alias.empty ())
62  produces <JetCollection>();
63  else
64  produces <JetCollection>().setBranchAlias(alias);
65  }
66 
67  template<class T>
69  const edm::EventSetup& fSetup)
70  {
71  // look for correctors
73  if (record.cacheIdentifier() != mCacheId)
74  { // need to renew cache
75  for (unsigned i = 0; i < mCorrectorNames.size(); i++)
76  {
78  record.get (mCorrectorNames [i], handle);
79  mCorrectors [i] = &*handle;
80  }
81  mCacheId = record.cacheIdentifier();
82  }
83  edm::Handle<JetCollection> jets; //Define Inputs
84  fEvent.getByToken (mInput, jets); //Get Inputs
85  std::auto_ptr<JetCollection> result (new JetCollection); //Corrected jets
86  typename JetCollection::const_iterator jet;
87  for (jet = jets->begin(); jet != jets->end(); jet++)
88  {
89  const T* referenceJet = &*jet;
90  int index = jet-jets->begin();
92  T correctedJet = *jet; //copy original jet
93  if (mVerbose)
94  std::cout<<"JetCorrectionProducer::produce-> original jet: "
95  <<jet->print()<<std::endl;
96  for (unsigned i = 0; i < mCorrectors.size(); ++i)
97  {
98  if ( !(mCorrectors[i]->vectorialCorrection()) ) {
99  // Scalar correction
100  double scale = 1.;
101  if (!(mCorrectors[i]->refRequired()))
102  scale = mCorrectors[i]->correction (*referenceJet,fEvent,fSetup);
103  else
104  scale = mCorrectors[i]->correction (*referenceJet,jetRef,fEvent,fSetup);
105  if (mVerbose)
106  std::cout<<"JetCorrectionProducer::produce-> Corrector # "
107  <<i<<", correction factor: "<<scale<<std::endl;
108  correctedJet.scaleEnergy (scale); // apply scalar correction
109  referenceJet = &correctedJet;
110  } else {
111  // Vectorial correction
112  JetCorrector::LorentzVector corrected;
113  double scale = mCorrectors[i]->correction (*referenceJet, jetRef,
114  fEvent, fSetup, corrected);
115  if (mVerbose)
116  std::cout<<"JetCorrectionProducer::produce-> Corrector # "
117  <<i<<", correction factor: "<<scale<<std::endl;
118  correctedJet.setP4( corrected ); // apply vectorial correction
119  referenceJet = &correctedJet;
120  }
121  }
122  if (mVerbose)
123  std::cout<<"JetCorrectionProducer::produce-> corrected jet: "
124  <<correctedJet.print ()<<std::endl;
125  result->push_back (correctedJet);
126  }
127  NumericSafeGreaterByPt<T> compJets;
128  // reorder corrected jets
129  std::sort (result->begin (), result->end (), compJets);
130  // put corrected jet collection into event
131  fEvent.put(result);
132  }
133 
134 }
135 
136 #endif
unsigned long long cacheIdentifier() const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
virtual void produce(edm::Event &, const edm::EventSetup &)
JetCorrectorParameters::Record record
Definition: classes.h:7
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
std::vector< std::string > mCorrectorNames
edm::EDGetTokenT< JetCollection > mInput
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
vector< PseudoJet > jets
tuple result
Definition: query.py:137
void get(HolderT &iHolder) const
tuple handle
Definition: patZpeak.py:22
JetCorrectionProducer(const edm::ParameterSet &fParameters)
const T & get() const
Definition: EventSetup.h:56
mVerbose(fConfig.getUntrackedParameter< bool >("verbose", false))
tuple cout
Definition: gather_cfg.py:121
volatile std::atomic< bool > shutdown_flag false
std::vector< const JetCorrector * > mCorrectors
long double T
tuple size
Write out results.
reco::Particle::LorentzVector LorentzVector
Definition: JetCorrector.h:23