17 #include "boost/bind.hpp"
44 namespace reco {
namespace tau {
50 double factor =
sqrt(vec.energy()*vec.energy()-mass*mass)/vec.P();
52 vec.px()*factor,vec.py()*factor,vec.pz()*factor,vec.energy());
72 std::unique_ptr<RecoTauQualityCuts>
qcuts_;
101 formula = formula.ReplaceAll(
"pT",
"x");
102 std::unique_ptr<TFormula>
function(
new TFormula(functionName.data(), formula.Data()));
103 int numParameter =
function->GetNpar();
104 for (
int idxParameter = 0; idxParameter < numParameter; ++idxParameter ) {
105 std::string parameterName = Form(
"par%i", idxParameter);
107 function->SetParameter(idxParameter, parameter);
129 qcuts_pset.addParameter<
double>(
"maxTrackChi2", 1.e+9);
130 qcuts_pset.addParameter<
double>(
"maxTransverseImpactParameter", 1.e+9);
131 qcuts_pset.addParameter<
double>(
"maxDeltaZ", 1.e+9);
132 qcuts_pset.addParameter<
double>(
"minTrackVertexWeight", -1.);
133 qcuts_pset.addParameter<
unsigned>(
"minTrackPixelHits", 0);
134 qcuts_pset.addParameter<
unsigned>(
"minTrackHits", 0);
153 if ( combineStrips_ ) {
172 std::set<size_t>& candIdsCurrentStrip,
bool& isCandAdded)
const
175 edm::LogPrint(
"RecoTauPiZeroStripPlugin3") <<
"<RecoTauPiZeroStripPlugin3::addCandsToStrip>:" ;
177 size_t numCands = cands.size();
178 for (
size_t candId = 0; candId < numCands; ++candId ) {
179 if ( (!candFlags[candId]) && candIdsCurrentStrip.find(candId) == candIdsCurrentStrip.end() ) {
183 if ( fabs(strip.
eta() - cand->eta()) < etaAssociationDistance_value &&
186 edm::LogPrint(
"RecoTauPiZeroStripPlugin3") <<
"--> adding PFCand #" << candId <<
" (" << cand.
id() <<
":" << cand.
key() <<
"): Et = " << cand->et() <<
", eta = " << cand->eta() <<
", phi = " << cand->phi() ;
191 candIdsCurrentStrip.insert(candId);
199 void markCandsInStrip(std::vector<bool>& candFlags,
const std::set<size_t>& candIds)
201 for ( std::set<size_t>::const_iterator candId = candIds.begin();
202 candId != candIds.end(); ++candId ) {
203 candFlags[*candId] =
true;
210 else if ( cand.gsfTrackRef().isNonnull() )
return reco::TrackBaseRef(cand.gsfTrackRef());
218 edm::LogPrint(
"RecoTauPiZeroStripPlugin3") <<
"<RecoTauPiZeroStripPlugin3::operator()>:" ;
234 for ( PFCandPtrs::iterator cand = candsVector.begin();
235 cand != candsVector.end(); ++cand ) {
237 edm::LogPrint(
"RecoTauPiZeroStripPlugin3") <<
"PFGamma #" << idx <<
" (" << cand->id() <<
":" << cand->key() <<
"): Et = " << (*cand)->et() <<
", eta = " << (*cand)->eta() <<
", phi = " << (*cand)->phi() ;
241 edm::LogPrint(
"RecoTauPiZeroStripPlugin3") <<
"--> assigning seedCandId = " << seedCands.size() ;
244 edm::LogPrint(
"RecoTauPiZeroStripPlugin3") <<
"track: Pt = " << candTrack->
pt() <<
" eta = " << candTrack->
eta() <<
", phi = " << candTrack->
phi() <<
", charge = " << candTrack->
charge() ;
247 <<
" chi2 = " << candTrack->
normalizedChi2() <<
", dPt/Pt = " << (candTrack->
ptError()/candTrack->
pt()) <<
")" ;
249 edm::LogPrint(
"RecoTauPiZeroStripPlugin3") <<
"ECAL Et: calibrated = " << (*cand)->ecalEnergy()*
sin((*cand)->theta()) <<
","
250 <<
" raw = " << (*cand)->rawEcalEnergy()*
sin((*cand)->theta()) ;
251 edm::LogPrint(
"RecoTauPiZeroStripPlugin3") <<
"HCAL Et: calibrated = " << (*cand)->hcalEnergy()*
sin((*cand)->theta()) <<
","
252 <<
" raw = " << (*cand)->rawHcalEnergy()*
sin((*cand)->theta()) ;
254 seedCands.push_back(*cand);
257 edm::LogPrint(
"RecoTauPiZeroStripPlugin3") <<
"--> assigning addCandId = " << addCands.size() ;
259 addCands.push_back(*cand);
264 std::vector<bool> seedCandFlags(seedCands.size());
265 std::vector<bool> addCandFlags(addCands.size());
267 std::set<size_t> seedCandIdsCurrentStrip;
268 std::set<size_t> addCandIdsCurrentStrip;
271 while ( idxSeed < seedCands.size() ) {
274 seedCandIdsCurrentStrip.clear();
275 addCandIdsCurrentStrip.clear();
278 strip->addDaughter(seedCands[idxSeed]);
279 seedCandIdsCurrentStrip.insert(idxSeed);
282 int stripBuildIteration = 0;
287 addCandsToStrip(*strip, seedCands, seedCandFlags, seedCandIdsCurrentStrip, isCandAdded);
289 addCandsToStrip(*strip, addCands, addCandFlags, addCandIdsCurrentStrip, isCandAdded);
293 ++stripBuildIteration;
297 if (
verbosity_ >= 2 )
edm::LogPrint(
"RecoTauPiZeroStripPlugin3") <<
"Building strip: Et = " << strip->et() <<
", eta = " << strip->eta() <<
", phi = " << strip->phi() ;
300 if ( strip->daughterPtr(0).isNonnull() ) strip->setVertex(strip->daughterPtr(0)->vertex());
301 output.push_back(strip);
304 markCandsInStrip(seedCandFlags, seedCandIdsCurrentStrip);
305 markCandsInStrip(addCandFlags, addCandIdsCurrentStrip);
307 if (
verbosity_ >= 2 )
edm::LogPrint(
"RecoTauPiZeroStripPlugin3") <<
"Discarding strip: Et = " << strip->et() <<
", eta = " << strip->eta() <<
", phi = " << strip->phi() ;
311 while ( idxSeed < seedCands.size() && seedCandFlags[idxSeed] ) {
320 output.sort(output.begin(), output.end(),
328 for ( PiZeroVector::const_iterator
first = output.begin();
330 for ( PiZeroVector::const_iterator
second =
first+1;
339 std::auto_ptr<RecoTauPiZero> combinedStrips(
346 for (
auto const& gamma :
first->daughterPtrVector()) {
347 combinedStrips->addDaughter(gamma);
349 for (
auto const& gamma :
second->daughterPtrVector()) {
350 combinedStrips->addDaughter(gamma);
353 if ( combinedStrips->daughterPtr(0).isNonnull() ) {
354 combinedStrips->setVertex(combinedStrips->daughterPtr(0)->vertex());
358 stripCombinations.push_back(combinedStrips);
363 output.transfer(output.end(), stripCombinations);
368 for ( PiZeroVector::iterator strip = output.begin();
369 strip != output.end(); ++strip ) {
370 double bendCorrEta = 0.;
371 double bendCorrPhi = 0.;
373 for (
auto const& gamma : strip->daughterPtrVector()) {
376 energySum += gamma->energy();
378 if ( energySum > 1.
e-2 ) {
383 strip->setBendCorrEta(bendCorrEta);
384 strip->setBendCorrPhi(bendCorrPhi);
387 return output.release();
T getParameter(std::string const &) const
InputIterator takeNElements(const InputIterator &begin, const InputIterator &end, size_t N)
virtual ~RecoTauPiZeroStripPlugin3()
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Sin< T >::type sin(const T &t)
bool updateStripAfterEachDaughter_
reco::VertexRef associatedVertex(const PFJet &jet) const
bool exists(std::string const ¶meterName) const
checks if a parameter exists
virtual double phi() const final
momentum azimuthal angle
double phi() const
azimuthal angle of momentum vector
void setEvent(const edm::Event &evt)
Load the vertices from the event.
return_type operator()(const reco::PFJet &) const override
Build a collection of piZeros from objects in the input jet.
bool isNonnull() const
Checks for non-null.
Jets made from PFObjects.
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
U second(std::pair< T, U > const &p)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
double eta() const
pseudorapidity of momentum vector
edm::RefToBase< reco::Track > TrackBaseRef
persistent reference to a Track, using views
double pt() const
track transverse momentum
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
void addCandsToStrip(RecoTauPiZero &, PFCandPtrs &, const std::vector< bool > &, std::set< size_t > &, bool &) const
bool applyElecTrackQcuts_
virtual void beginEvent() override
Hook called at the beginning of the event.
std::unique_ptr< const TFormula > etaAssociationDistance_
std::unique_ptr< RecoTauQualityCuts > qcuts_
std::vector< reco::PFCandidatePtr > PFCandPtrs
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
double deltaPhi(double phi1, double phi2)
static const TrackGhostTrackState * getTrack(const BasicGhostTrackState *basic)
double minGammaEtStripSeed_
const edm::Event * evt() const
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
ParameterSet const & getParameterSet(std::string const &) const
int maxStripBuildIterations_
tuple idx
DEBUGGING if hasattr(process,"trackMonIterativeTracking2012"): print "trackMonIterativeTracking2012 D...
double combinatoricStripMassHypo_
ProductID id() const
Accessor for product ID.
math::XYZTLorentzVector LorentzVector
Lorentz vector.
std::vector< int > inputPdgIds_
int numberOfValidTrackerHits() const
std::auto_ptr< PiZeroVector > return_type
RecoTauPiZeroStripPlugin3(const edm::ParameterSet &, edm::ConsumesCollector &&iC)
boost::ptr_vector< RecoTauPiZero > PiZeroVector
std::unique_ptr< const TFormula > phiAssociationDistance_
RecoTauVertexAssociator vertexAssociator_
AddFourMomenta p4Builder_
void addDaughter(const CandidatePtr &)
add a daughter via a reference
int numberOfValidPixelHits() const
double minGammaEtStripAdd_
math::XYZPoint Point
point in the space
int charge() const
track electric charge
#define DEFINE_EDM_PLUGIN(factory, type, name)
void set(reco::Candidate &c) const
set up a candidate
virtual double eta() const final
momentum pseudorapidity
double energySum(const DataFrame &df, int fs, int ls)
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
virtual double pt() const final
transverse momentum