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 
12 
13 #define init_param(type, varname) varname (ps.getParameter< type >( #varname ))
14 
15 namespace fftjetcms {
16 
18 {
19  return true;
20 }
21 
22 
24  : inputLabel(ps.getParameter<edm::InputTag>("src")),
25  init_param(std::string, outputLabel),
26  jetType(parseJetType(ps.getParameter<std::string>("jetType"))),
27  init_param(bool, doPVCorrection),
28  init_param(edm::InputTag, srcPVs),
29  etaDependentMagnutideFactors(
30  ps.getParameter<std::vector<double> >(
31  "etaDependentMagnutideFactors")),
32  init_param(edm::ParameterSet, anomalous),
33  init_param(bool, insertCompleteEvent),
34  init_param(double, completeEventScale),
35  vertex_(0.0, 0.0, 0.0)
36 {
38  throw cms::Exception("FFTJetBadConfig")
39  << "Bad scale for the complete event : must be positive"
40  << std::endl;
41 }
42 
43 
45 {
47 }
48 
49 
51 {
52  // Figure out if we are going to perform the vertex adjustment
53  const bool adjustForVertex = doPVCorrection && jetType == CALOJET;
54 
55  // Figure out the vertex
56  if (adjustForVertex)
57  {
59  iEvent.getByLabel(srcPVs, pvCollection);
60  if (pvCollection->empty())
61  vertex_ = reco::Particle::Point(0.0, 0.0, 0.0);
62  else
63  vertex_ = pvCollection->begin()->position();
64  }
65 
66  // Get the input collection
68 
69  // Create the set of 4-vectors needed by the algorithm
70  eventData.clear();
71  candidateIndex.clear();
72  unsigned index = 0;
75  it != end; ++it, ++index)
76  {
77  const reco::Candidate& item(*it);
78  if (anomalous(item))
79  continue;
80  if (edm::isNotFinite(item.pt()))
81  continue;
82 
83  if (adjustForVertex)
84  {
85  const CaloTower& tower(dynamic_cast<const CaloTower&>(item));
86  eventData.push_back(VectorLike(tower.p4(vertex_)));
87  }
88  else
89  eventData.push_back(item.p4());
90  candidateIndex.push_back(index);
91  }
92  assert(eventData.size() == candidateIndex.size());
93 }
94 
95 
97 {
98  // It is a bug to call this function before defining the energy flow grid
99  assert(energyFlow.get());
100 
101  fftjet::Grid2d<Real>& g(*energyFlow);
102  g.reset();
103 
104  const unsigned nInputs = eventData.size();
105  const VectorLike* inp = nInputs ? &eventData[0] : 0;
106  for (unsigned i=0; i<nInputs; ++i)
107  {
108  const VectorLike& item(inp[i]);
109  g.fill(item.Eta(), item.Phi(), item.Et());
110  }
111 
112  if (!etaDependentMagnutideFactors.empty())
113  {
114  if (etaDependentMagnutideFactors.size() != g.nEta())
115  throw cms::Exception("FFTJetBadConfig")
116  << "ERROR in FFTJetInterface::discretizeEnergyFlow() :"
117  " number of elements in the \"etaDependentMagnutideFactors\""
118  " vector is inconsistent with the grid binning"
119  << std::endl;
120  g.scaleData(&etaDependentMagnutideFactors[0],
122  }
123 }
124 
125 }
int i
Definition: DBlmapReader.cc:9
std::auto_ptr< fftjet::Grid2d< fftjetcms::Real > > energyFlow
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
const edm::InputTag inputLabel
edm::Handle< reco::CandidateView > inputCollection
math::PtEtaPhiMLorentzVector p4(double vtxZ) const
Definition: CaloTower.cc:105
void loadInputCollection(const edm::Event &)
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:243
std::vector< unsigned > candidateIndex
bool isNotFinite(T x)
Definition: isFinite.h:10
reco::Particle::Point vertex_
math::XYZTLorentzVector VectorLike
const std::vector< double > etaDependentMagnutideFactors
math::XYZPoint Point
point in the space
Definition: Particle.h:29
#define end
Definition: vmac.h:38
const AnomalousTower anomalous
JetType parseJetType(const std::string &name)
Definition: JetType.cc:5
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
bool storeInSinglePrecision() const
#define init_param(type, varname)
std::vector< fftjetcms::VectorLike > eventData
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector