CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FFTJetInterface.cc
Go to the documentation of this file.
1 #include <cassert>
2 
4 
7 
10 
11 #define init_param(type, varname) varname (ps.getParameter< type >( #varname ))
12 
13 namespace fftjetcms {
14 
16 {
17  return true;
18 }
19 
20 
22  : edm::EDProducer(),
23  inputLabel(ps.getParameter<edm::InputTag>("src")),
24  init_param(std::string, outputLabel),
25  jetType(parseJetType(ps.getParameter<std::string>("jetType"))),
26  init_param(bool, doPVCorrection),
27  init_param(edm::InputTag, srcPVs),
28  etaDependentMagnutideFactors(
29  ps.getParameter<std::vector<double> >(
30  "etaDependentMagnutideFactors")),
31  init_param(edm::ParameterSet, anomalous),
32  init_param(bool, insertCompleteEvent),
33  init_param(double, completeEventScale),
34  vertex_(0.0, 0.0, 0.0)
35 {
37  throw cms::Exception("FFTJetBadConfig")
38  << "Bad scale for the complete event : must be positive"
39  << std::endl;
40 
41  inputToken = consumes<reco::CandidateView>(inputLabel);
42 
43  if (doPVCorrection && jetType == CALOJET)
44  srcPVsToken = consumes<reco::VertexCollection>(srcPVs);
45 }
46 
47 
49 {
51 }
52 
53 
55 {
56  // Figure out if we are going to perform the vertex adjustment
57  const bool adjustForVertex = doPVCorrection && jetType == CALOJET;
58 
59  // Figure out the vertex
60  if (adjustForVertex)
61  {
63  iEvent.getByToken(srcPVsToken, pvCollection);
64  if (pvCollection->empty())
65  vertex_ = reco::Particle::Point(0.0, 0.0, 0.0);
66  else
67  vertex_ = pvCollection->begin()->position();
68  }
69 
70  // Get the input collection
72 
73  // Create the set of 4-vectors needed by the algorithm
74  eventData.clear();
75  candidateIndex.clear();
76  unsigned index = 0;
79  it != end; ++it, ++index)
80  {
81  const reco::Candidate& item(*it);
82  if (anomalous(item))
83  continue;
84  if (edm::isNotFinite(item.pt()))
85  continue;
86 
87  if (adjustForVertex)
88  {
89  const CaloTower& tower(dynamic_cast<const CaloTower&>(item));
90  eventData.push_back(VectorLike(tower.p4(vertex_)));
91  }
92  else
93  eventData.push_back(item.p4());
94  candidateIndex.push_back(index);
95  }
96  assert(eventData.size() == candidateIndex.size());
97 }
98 
99 
101 {
102  // It is a bug to call this function before defining the energy flow grid
103  assert(energyFlow.get());
104 
105  fftjet::Grid2d<Real>& g(*energyFlow);
106  g.reset();
107 
108  const unsigned nInputs = eventData.size();
109  const VectorLike* inp = nInputs ? &eventData[0] : 0;
110  for (unsigned i=0; i<nInputs; ++i)
111  {
112  const VectorLike& item(inp[i]);
113  g.fill(item.Eta(), item.Phi(), item.Et());
114  }
115 
116  if (!etaDependentMagnutideFactors.empty())
117  {
118  if (etaDependentMagnutideFactors.size() != g.nEta())
119  throw cms::Exception("FFTJetBadConfig")
120  << "ERROR in FFTJetInterface::discretizeEnergyFlow() :"
121  " number of elements in the \"etaDependentMagnutideFactors\""
122  " vector is inconsistent with the grid binning"
123  << std::endl;
124  g.scaleData(&etaDependentMagnutideFactors[0],
126  }
127 }
128 
129 }
edm::EDGetTokenT< reco::CandidateView > inputToken
int i
Definition: DBlmapReader.cc:9
std::auto_ptr< fftjet::Grid2d< fftjetcms::Real > > energyFlow
const edm::InputTag inputLabel
edm::Handle< reco::CandidateView > inputCollection
math::PtEtaPhiMLorentzVector p4(double vtxZ) const
Definition: CaloTower.cc:129
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
void loadInputCollection(const edm::Event &)
tuple jetType
JET
Definition: autophobj.py:147
const edm::InputTag srcPVs
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
virtual float pt() const =0
transverse momentum
int iEvent
Definition: GenABIO.cc:230
std::vector< unsigned > candidateIndex
bool isNotFinite(T x)
Definition: isFinite.h:10
reco::Particle::Point vertex_
math::XYZTLorentzVector VectorLike
const std::vector< double > etaDependentMagnutideFactors
edm::EDGetTokenT< reco::VertexCollection > srcPVsToken
math::XYZPoint Point
point in the space
Definition: Particle.h:31
#define end
Definition: vmac.h:37
const AnomalousTower anomalous
JetType parseJetType(const std::string &name)
Definition: JetType.cc:5
bool storeInSinglePrecision() const
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
#define init_param(type, varname)
std::vector< fftjetcms::VectorLike > eventData
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector