CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
OniaPhotonConversionProducer Class Reference

#include <OniaPhotonConversionProducer.h>

Inheritance diagram for OniaPhotonConversionProducer:
edm::stream::EDProducer<> edm::stream::EDProducerBase edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 OniaPhotonConversionProducer (const edm::ParameterSet &ps)
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
- Public Member Functions inherited from edm::stream::EDProducerBase
 EDProducerBase ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDProducerBase ()
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector< edm::ProductResolverIndex > const & indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription const &)> registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, std::unordered_multimap< std::string, edm::ProductResolverIndex > const &iIndicies, std::string const &moduleLabel)
 
virtual ~ProducerBase () noexcept(false)
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

bool CheckPi0 (const reco::Conversion &, const reco::PFCandidateCollection &, bool &)
 
bool checkTkVtxCompatibility (const reco::Conversion &, const reco::VertexCollection &)
 
reco::Candidate::LorentzVector convertVector (const math::XYZTLorentzVectorF &)
 
virtual void endStream ()
 
bool foundCompatibleInnerHits (const reco::HitPattern &hitPatA, const reco::HitPattern &hitPatB)
 
bool HighpuritySubset (const reco::Conversion &, const reco::VertexCollection &)
 
pat::CompositeCandidatemakePhotonCandidate (const reco::Conversion &)
 
int PackFlags (const reco::Conversion &, bool, bool, bool, bool, bool)
 
virtual void produce (edm::Event &event, const edm::EventSetup &esetup)
 
void removeDuplicates (reco::ConversionCollection &)
 
const reco::PFCandidateCollection selectPFPhotons (const reco::PFCandidateCollection &)
 

Private Attributes

double _minDistanceOfApproachMaxCut
 
double _minDistanceOfApproachMinCut
 
double _trackchi2Cut
 
double _vertexChi2ProbCut
 
int convAlgo_
 
edm::EDGetTokenT< reco::ConversionCollectionconvCollectionToken_
 
std::vector< int > convQuality_
 
std::string convSelectionCuts_
 
edm::EDGetTokenT< reco::PFCandidateCollectionpfCandidateCollectionToken_
 
std::vector< double > pi0LargeWindow_
 
bool pi0OnlineSwitch_
 
std::vector< double > pi0SmallWindow_
 
uint32_t sigmaTkVtxComp_
 
edm::EDGetTokenT< reco::VertexCollectionthePVsToken_
 
uint32_t TkMinNumOfDOF_
 
bool wantCompatibleInnerHits_
 
bool wantHighpurity_
 
bool wantTkVtxCompatibility_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
typedef CacheContexts< T... > CacheTypes
 
typedef CacheTypes::GlobalCache GlobalCache
 
typedef AbilityChecker< T... > HasAbility
 
typedef CacheTypes::LuminosityBlockCache LuminosityBlockCache
 
typedef LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCacheLuminosityBlockContext
 
typedef CacheTypes::LuminosityBlockSummaryCache LuminosityBlockSummaryCache
 
typedef CacheTypes::RunCache RunCache
 
typedef RunContextT< RunCache, GlobalCacheRunContext
 
typedef CacheTypes::RunSummaryCache RunSummaryCache
 
- Public Types inherited from edm::stream::EDProducerBase
typedef EDProducerAdaptorBase ModuleType
 
- Public Types inherited from edm::ProducerBase
typedef ProductRegistryHelper::TypeLabelList TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::stream::EDProducerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Select photon conversions and produce a conversion candidate collection

Definition at line 36 of file OniaPhotonConversionProducer.h.

Constructor & Destructor Documentation

OniaPhotonConversionProducer::OniaPhotonConversionProducer ( const edm::ParameterSet ps)
explicit

Definition at line 55 of file OniaPhotonConversionProducer.cc.

References patPFMETCorrections_cff::algo, edm::ParameterSet::getParameter(), and AlCaHLTBitMon_QueryRunRegistry::string.

55  {
56  convCollectionToken_ = consumes<reco::ConversionCollection>(ps.getParameter<edm::InputTag>("conversions"));
57  thePVsToken_ = consumes<reco::VertexCollection>(ps.getParameter<edm::InputTag>("primaryVertexTag"));
58 
59  wantTkVtxCompatibility_ = ps.getParameter<bool>("wantTkVtxCompatibility");
60  sigmaTkVtxComp_ = ps.getParameter<uint32_t>("sigmaTkVtxComp");
61  wantCompatibleInnerHits_ = ps.getParameter<bool>("wantCompatibleInnerHits");
62  TkMinNumOfDOF_ = ps.getParameter<uint32_t>("TkMinNumOfDOF");
63 
64  wantHighpurity_ = ps.getParameter<bool>("wantHighpurity");
65  _vertexChi2ProbCut = ps.getParameter<double>("vertexChi2ProbCut");
66  _trackchi2Cut = ps.getParameter<double>("trackchi2Cut");
67  _minDistanceOfApproachMinCut = ps.getParameter<double>("minDistanceOfApproachMinCut");
68  _minDistanceOfApproachMaxCut = ps.getParameter<double>("minDistanceOfApproachMaxCut");
69 
70  pfCandidateCollectionToken_ = consumes<reco::PFCandidateCollection>(ps.getParameter<edm::InputTag>("pfcandidates"));
71 
72  pi0OnlineSwitch_ = ps.getParameter<bool>("pi0OnlineSwitch");
73  pi0SmallWindow_ = ps.getParameter<std::vector<double> >("pi0SmallWindow");
74  pi0LargeWindow_ = ps.getParameter<std::vector<double> >("pi0LargeWindow");
75 
76  std::string algo = ps.getParameter<std::string>("convAlgo");
77  convAlgo_ = StringToEnumValue<reco::Conversion::ConversionAlgorithm>(algo);
78 
79  std::vector<std::string> qual = ps.getParameter<std::vector<std::string> >("convQuality");
80  if( qual[0] != "" ) convQuality_ =StringToEnumValue<reco::Conversion::ConversionQuality>(qual);
81 
82  convSelectionCuts_ = ps.getParameter<std::string>("convSelection");
83  produces<pat::CompositeCandidateCollection>("conversions");
84 }
T getParameter(std::string const &) const
edm::EDGetTokenT< reco::VertexCollection > thePVsToken_
edm::EDGetTokenT< reco::ConversionCollection > convCollectionToken_
edm::EDGetTokenT< reco::PFCandidateCollection > pfCandidateCollectionToken_

Member Function Documentation

bool OniaPhotonConversionProducer::CheckPi0 ( const reco::Conversion conv,
const reco::PFCandidateCollection photons,
bool &  pizero_rejected 
)
private

Definition at line 328 of file OniaPhotonConversionProducer.cc.

References muons2muons_cfi::photon, and reco::Conversion::refittedPair4Momentum().

329  {
330  // 2 windows are defined for Pi0 rejection, Conversions that, paired with others photons from the event, have an
331  // invariant mass inside the "small" window will be pizero_rejected and those that falls in the large window will
332  // be CheckPi0.
333  bool check_small = false;
334  bool check_large = false;
335 
336  float small1 = pi0SmallWindow_[0];
337  float small2 = pi0SmallWindow_[1];
338  float large1 = pi0LargeWindow_[0];
339  float large2 = pi0LargeWindow_[1];
340  for (reco::PFCandidateCollection::const_iterator photon = photons.begin(); photon!=photons.end(); ++photon) {
341  float inv = (conv.refittedPair4Momentum() + photon->p4()).M();
342  if (inv > large1 && inv < large2) {
343  check_large = true;
344  if (inv > small1 && inv < small2) {
345  check_small = true;
346  break;
347  }
348  }
349  }
350  pizero_rejected = check_small;
351  return check_large;
352 }
math::XYZTLorentzVectorF refittedPair4Momentum() const
Conversion track pair 4-momentum from the tracks refitted with vertex constraint. ...
Definition: Conversion.cc:235
bool OniaPhotonConversionProducer::checkTkVtxCompatibility ( const reco::Conversion conv,
const reco::VertexCollection priVtxs 
)
private

Definition at line 233 of file OniaPhotonConversionProducer.cc.

References begin, KineDebug3::count(), reco::Vertex::covariance(), reco::TrackBase::dz(), reco::TrackBase::dzError(), end, training_settings::idx, lt_(), reco::Vertex::position(), edm::second(), findQualityFiles::size, mathSSE::sqrt(), reco::Conversion::tracks(), and badGlobalMuonTaggersAOD_cff::vtx.

233  {
234  std::vector< std::pair< double, short> > idx[2];
235  short ik=-1;
236  BOOST_FOREACH(edm::RefToBase<reco::Track> tk, conv.tracks()){
237  ik++;
238  short count=-1;
239  BOOST_FOREACH(const reco::Vertex& vtx,priVtxs){
240  count++;
241 
242  double dz_= tk->dz(vtx.position());
243  double dzError_=tk->dzError();
244  dzError_=sqrt(dzError_*dzError_+vtx.covariance(2,2));
245 
246  if(fabs(dz_)/dzError_ > sigmaTkVtxComp_) continue;
247 
248  idx[ik].push_back(std::pair<double,short>(fabs(dz_),count));
249  }
250  if(idx[ik].size()==0) {return false;}
251 
252  std::stable_sort(idx[ik].begin(),idx[ik].end(),lt_);
253  }
254  if (idx[0][0].second==idx[1][0].second || idx[0][1].second==idx[1][0].second || idx[0][0].second==idx[1][1].second) return true;
255  return false;
256 }
size
Write out results.
double covariance(int i, int j) const
(i, j)-th element of error matrix, i, j = 0, ... 2
Definition: Vertex.h:130
const Point & position() const
position
Definition: Vertex.h:109
U second(std::pair< T, U > const &p)
T sqrt(T t)
Definition: SSEVec.h:18
#define end
Definition: vmac.h:37
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...
Definition: TrackBase.h:604
double dzError() const
error on dz
Definition: TrackBase.h:809
bool lt_(std::pair< double, short > a, std::pair< double, short > b)
#define begin
Definition: vmac.h:30
std::vector< edm::RefToBase< reco::Track > > const & tracks() const
vector of track to base references
Definition: Conversion.cc:176
reco::Candidate::LorentzVector OniaPhotonConversionProducer::convertVector ( const math::XYZTLorentzVectorF v)
private

Definition at line 354 of file OniaPhotonConversionProducer.cc.

354  {
355  return reco::Candidate::LorentzVector(v.x(),v.y(), v.z(), v.t());
356 }
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
void OniaPhotonConversionProducer::endStream ( )
privatevirtual

Reimplemented from edm::stream::EDProducerBase.

Definition at line 359 of file OniaPhotonConversionProducer.cc.

References DEFINE_FWK_MODULE.

359  {
360 }
bool OniaPhotonConversionProducer::foundCompatibleInnerHits ( const reco::HitPattern hitPatA,
const reco::HitPattern hitPatB 
)
private

Definition at line 258 of file OniaPhotonConversionProducer.cc.

References KineDebug3::count(), reco::HitPattern::getHitPattern(), reco::HitPattern::getLayer(), reco::HitPattern::getSubStructure(), reco::HitPattern::getTrackerMonoStereo(), HighpuritySubset(), mps_fire::i, reco::HitPattern::numberOfHits(), reco::HitPattern::trackerHitFilter(), and reco::HitPattern::validHitFilter().

258  {
259  size_t count=0;
260  uint32_t oldSubStr=0;
261  for (int i=0; i<hitPatA.numberOfHits(reco::HitPattern::HitCategory::TRACK_HITS) && count<2; i++) {
262  uint32_t hitA = hitPatA.getHitPattern(reco::HitPattern::HitCategory::TRACK_HITS,i);
263  if (!hitPatA.validHitFilter(hitA) || !hitPatA.trackerHitFilter(hitA)) continue;
264 
265  if(hitPatA.getSubStructure(hitA)==oldSubStr && hitPatA.getLayer(hitA)==oldSubStr)
266  continue;
267 
268  if(hitPatB.getTrackerMonoStereo(reco::HitPattern::HitCategory::TRACK_HITS,hitPatA.getSubStructure(hitA),hitPatA.getLayer(hitA)) != 0)
269  return true;
270 
271  oldSubStr=hitPatA.getSubStructure(hitA);
272  count++;
273  }
274  return false;
275 }
static uint32_t getLayer(uint16_t pattern)
Definition: HitPattern.h:700
uint16_t getTrackerMonoStereo(HitCategory category, uint16_t substr, uint16_t layer) const
Definition: HitPattern.cc:454
static bool validHitFilter(uint16_t pattern)
Definition: HitPattern.h:787
static uint32_t getSubStructure(uint16_t pattern)
Definition: HitPattern.h:691
static bool trackerHitFilter(uint16_t pattern)
Definition: HitPattern.h:677
uint16_t getHitPattern(HitCategory category, int position) const
Definition: HitPattern.h:515
int numberOfHits(HitCategory category) const
Definition: HitPattern.h:807
bool OniaPhotonConversionProducer::HighpuritySubset ( const reco::Conversion conv,
const reco::VertexCollection priVtxs 
)
private

Definition at line 278 of file OniaPhotonConversionProducer.cc.

References reco::TrackBase::charge(), reco::Vertex::chi2(), ChiSquaredProbability(), reco::Conversion::conversionVertex(), reco::Conversion::distOfMinimumApproach(), reco::TrackBase::dxy(), reco::TrackBase::dxyError(), mps_fire::i, reco::Vertex::ndof(), reco::TrackBase::ndof(), reco::TrackBase::normalizedChi2(), reco::Vertex::position(), reco::Conversion::tracks(), badGlobalMuonTaggersAOD_cff::vtx, and reco::Conversion::zOfPrimaryVertexFromTracks().

Referenced by foundCompatibleInnerHits().

278  {
279  // select high purity conversions our way:
280  // vertex chi2 cut
282 
283  // d0 cut
284  // Find closest primary vertex
285  int closest_pv_index = 0;
286  int i=0;
287  BOOST_FOREACH(const reco::Vertex& vtx,priVtxs){
288  if( conv.zOfPrimaryVertexFromTracks( vtx.position() ) < conv.zOfPrimaryVertexFromTracks( priVtxs[closest_pv_index].position() ) ) closest_pv_index = i;
289  i++;
290  }
291  // Now check impact parameter wtr with the just found closest primary vertex
292  BOOST_FOREACH(const edm::RefToBase<reco::Track> tk, conv.tracks()) if(-tk->dxy(priVtxs[closest_pv_index].position())*tk->charge()/tk->dxyError()<0) return false;
293 
294  // chi2 of single tracks
295  BOOST_FOREACH(const edm::RefToBase<reco::Track> tk, conv.tracks()) if(tk->normalizedChi2() > _trackchi2Cut) return false;
296 
297  // dof for each track
298  BOOST_FOREACH(const edm::RefToBase<reco::Track> tk, conv.tracks()) if(tk->ndof()< TkMinNumOfDOF_) return false;
299 
300  // distance of approach cut
302 
303  return true;
304 }
const reco::Vertex & conversionVertex() const
returns the reco conversion vertex
Definition: Conversion.h:97
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Definition: TrackBase.h:556
double dxyError() const
error on dxy
Definition: TrackBase.h:791
double zOfPrimaryVertexFromTracks(const math::XYZPoint &myBeamSpot=math::XYZPoint()) const
Definition: Conversion.h:145
double distOfMinimumApproach() const
Definition: Conversion.h:125
const Point & position() const
position
Definition: Vertex.h:109
double ndof() const
number of degrees of freedom of the fit
Definition: TrackBase.h:550
double chi2() const
chi-squares
Definition: Vertex.h:98
float ChiSquaredProbability(double chiSquared, double nrDOF)
double ndof() const
Definition: Vertex.h:105
int charge() const
track electric charge
Definition: TrackBase.h:562
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:586
std::vector< edm::RefToBase< reco::Track > > const & tracks() const
vector of track to base references
Definition: Conversion.cc:176
pat::CompositeCandidate * OniaPhotonConversionProducer::makePhotonCandidate ( const reco::Conversion conv)
private

Definition at line 306 of file OniaPhotonConversionProducer.cc.

References pat::PATObject< ObjectType >::addUserData(), reco::Conversion::conversionVertex(), reco::Vertex::position(), reco::Conversion::refittedPair4Momentum(), reco::LeafCandidate::setP4(), reco::LeafCandidate::setVertex(), and reco::Conversion::tracks().

306  {
307 
309  photonCand->setP4(convertVector(conv.refittedPair4Momentum()));
310  photonCand->setVertex(conv.conversionVertex().position());
311 
312  photonCand->addUserData<reco::Track>( "track0", *conv.tracks()[0] );
313  photonCand->addUserData<reco::Track>( "track1", *conv.tracks()[1] );
314 
315  return photonCand;
316 
317 }
const reco::Vertex & conversionVertex() const
returns the reco conversion vertex
Definition: Conversion.h:97
Analysis-level particle class.
reco::Candidate::LorentzVector convertVector(const math::XYZTLorentzVectorF &)
const Point & position() const
position
Definition: Vertex.h:109
virtual void setVertex(const Point &vertex)
set vertex
math::XYZTLorentzVectorF refittedPair4Momentum() const
Conversion track pair 4-momentum from the tracks refitted with vertex constraint. ...
Definition: Conversion.cc:235
virtual void setP4(const LorentzVector &p4) final
set 4-momentum
void addUserData(const std::string &label, const T &data, bool transientOnly=false, bool overwrite=false)
Definition: PATObject.h:309
std::vector< edm::RefToBase< reco::Track > > const & tracks() const
vector of track to base references
Definition: Conversion.cc:176
int OniaPhotonConversionProducer::PackFlags ( const reco::Conversion conv,
bool  flagTkVtxCompatibility,
bool  flagCompatibleInnerHits,
bool  flagHighpurity,
bool  pizero_rejected,
bool  large_pizero_window 
)
private

Definition at line 182 of file OniaPhotonConversionProducer.cc.

References reco::Conversion::algo(), flags, lumiQueryAPI::q, and reco::Conversion::quality().

184  {
185  int flags = 0;
186  if ( flagTkVtxCompatibility ) flags += 1;
187  if ( flagCompatibleInnerHits ) flags += 2;
188  if ( flagHighpurity ) flags += 4;
189  if ( pizero_rejected ) flags += 8;
190  if ( large_pizero_window ) flags += 16;
191 
192  flags += (conv.algo()*32);
193  int q_mask = 0;
194  std::vector<std::string> s_quals;
195  s_quals.push_back("generalTracksOnly");
196  s_quals.push_back("arbitratedEcalSeeded");
197  s_quals.push_back("arbitratedMerged");
198  s_quals.push_back("arbitratedMergedEcalGeneral");
199  s_quals.push_back("highPurity");
200  s_quals.push_back("highEfficiency");
201  s_quals.push_back("ecalMatched1Track");
202  s_quals.push_back("ecalMatched2Track");
203  std::vector<int> i_quals = StringToEnumValue<reco::Conversion::ConversionQuality>(s_quals);
204  for (std::vector<int>::const_iterator qq = i_quals.begin(); qq!=i_quals.end(); ++qq) {
206  if (conv.quality(q)) q_mask = *qq;
207  }
208  flags += (q_mask*32*8);
209  return flags;
210 }
ConversionAlgorithm algo() const
Definition: Conversion.h:223
bool quality(ConversionQuality q) const
Definition: Conversion.h:181
std::vector< Variable::Flags > flags
Definition: MVATrainer.cc:135
void OniaPhotonConversionProducer::produce ( edm::Event event,
const edm::EventSetup esetup 
)
privatevirtual

Definition at line 86 of file OniaPhotonConversionProducer.cc.

References pat::PATObject< ObjectType >::addUserInt(), conv, RemoveAddSevLevel::flag, flags, eostools::move(), edm::Handle< T >::product(), and lumiQueryAPI::q.

86  {
87 
88  std::unique_ptr<reco::ConversionCollection> outCollection(new reco::ConversionCollection);
89  std::unique_ptr<pat::CompositeCandidateCollection> patoutCollection(new pat::CompositeCandidateCollection);
90 
92  event.getByToken(thePVsToken_, priVtxs);
93 
95  event.getByToken(convCollectionToken_,pConv);
96 
98  event.getByToken(pfCandidateCollectionToken_,pfcandidates);
99 
100  const reco::PFCandidateCollection pfphotons = selectPFPhotons(*pfcandidates);
101 
103 
104  for(reco::ConversionCollection::const_iterator conv = pConv->begin(); conv != pConv->end(); ++conv){
105 
106  if (! ( *convSelection_)(*conv)){
107  continue; // selection string
108  }
109  if (convAlgo_ != 0 && conv->algo()!= convAlgo_){
110  continue; // select algorithm
111  }
112  if(convQuality_.size() > 0){
113  bool flagsok=true;
114  for (std::vector<int>::const_iterator flag = convQuality_.begin(); flag!=convQuality_.end(); ++flag){
116  if (!conv->quality(q)) {
117  flagsok=false;
118  break;
119  }
120  }
121  if (!flagsok){
122  continue;
123  }
124  }
125  outCollection->push_back(*conv);
126  }
127 
128  removeDuplicates(*outCollection);
129 
130  for (reco::ConversionCollection::const_iterator conv = outCollection->begin(); conv != outCollection->end(); ++conv){
131 
132  bool flag1 = true;
133  bool flag2 = true;
134  bool flag3 = true;
135  bool flag4 = true;
136 
137  // The logic implies that by default this flags are true and if the check are not wanted conversions are saved.
138  // If checks are required and failed then don't save the conversion.
139 
140  bool flagTkVtxCompatibility = true;
141  if (! checkTkVtxCompatibility(*conv,*priVtxs.product())) {
142  flagTkVtxCompatibility = false;
144  flag1 = false;
145  }
146  }
147  bool flagCompatibleInnerHits = false;
148  if (conv->tracks().size()==2) {
149  reco::HitPattern hitPatA=conv->tracks().at(0)->hitPattern();
150  reco::HitPattern hitPatB=conv->tracks().at(1)->hitPattern();
151  if ( foundCompatibleInnerHits(hitPatA,hitPatB) && foundCompatibleInnerHits(hitPatB,hitPatA) ) flagCompatibleInnerHits = true;
152  }
153  if (wantCompatibleInnerHits_ && ! flagCompatibleInnerHits) {
154  flag2 = false;
155  }
156  bool flagHighpurity = true;
157  if (!HighpuritySubset(*conv,*priVtxs.product())) {
158  flagHighpurity = false;
159  if (wantHighpurity_) {
160  flag3 = false;
161  }
162  }
163  bool pizero_rejected = false;
164  bool large_pizero_window = CheckPi0(*conv, pfphotons, pizero_rejected);
165  if (pi0OnlineSwitch_ && pizero_rejected) {
166  flag4 = false;
167  }
168 
169  int flags = 0;
170  if (flag1 && flag2 && flag3 && flag4){
171  flags = PackFlags(*conv,flagTkVtxCompatibility,flagCompatibleInnerHits,flagHighpurity,pizero_rejected,large_pizero_window);
173  pat_conv->addUserInt("flags",flags);
174  patoutCollection->push_back(*pat_conv);
175  }
176  }
177  event.put(std::move(patoutCollection),"conversions");
178 
179  delete convSelection_;
180 }
Analysis-level particle class.
bool foundCompatibleInnerHits(const reco::HitPattern &hitPatA, const reco::HitPattern &hitPatB)
static HepMC::IO_HEPEVT conv
bool HighpuritySubset(const reco::Conversion &, const reco::VertexCollection &)
std::vector< Variable::Flags > flags
Definition: MVATrainer.cc:135
std::vector< Conversion > ConversionCollection
collectin of Conversion objects
Definition: ConversionFwd.h:9
bool checkTkVtxCompatibility(const reco::Conversion &, const reco::VertexCollection &)
int PackFlags(const reco::Conversion &, bool, bool, bool, bool, bool)
edm::EDGetTokenT< reco::VertexCollection > thePVsToken_
void addUserInt(const std::string &label, int32_t data, const bool overwrite=false)
Set user-defined int.
Definition: PATObject.h:855
const reco::PFCandidateCollection selectPFPhotons(const reco::PFCandidateCollection &)
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
edm::EDGetTokenT< reco::ConversionCollection > convCollectionToken_
std::vector< CompositeCandidate > CompositeCandidateCollection
T const * product() const
Definition: Handle.h:81
edm::EDGetTokenT< reco::PFCandidateCollection > pfCandidateCollectionToken_
bool CheckPi0(const reco::Conversion &, const reco::PFCandidateCollection &, bool &)
pat::CompositeCandidate * makePhotonCandidate(const reco::Conversion &)
void removeDuplicates(reco::ConversionCollection &)
def move(src, dest)
Definition: eostools.py:510
void OniaPhotonConversionProducer::removeDuplicates ( reco::ConversionCollection c)
private

Put in out collection only those conversion candidates that are not sharing tracks. If sharing, keep the one with the best chi2.

Definition at line 215 of file OniaPhotonConversionProducer.cc.

References ConversionEqualByTrack(), and ConversionLessByChi2().

215  {
216  // first sort from high to low chi2 prob
217  std::sort(c.begin(),c.end(),ConversionLessByChi2);
218  int iter1 = 0;
219  // Cycle over all the elements of the collection and compare to all the following,
220  // if two elements have at least one track in common delete the element with the lower chi2
221  while(iter1 < (( (int) c.size() ) - 1) ){
222  int iter2 = iter1+1;
223  while( iter2 < (int) c.size()) if(ConversionEqualByTrack( c[iter1], c[iter2] ) ){
224  c.erase( c.begin() + iter2 );
225  }else{
226  iter2++; // Increment index only if this element is no duplicate.
227  // If it is duplicate check again the same index since the vector rearranged elements index after erasing
228  }
229  iter1++;
230  }
231 }
bool ConversionLessByChi2(const reco::Conversion &c1, const reco::Conversion &c2)
bool ConversionEqualByTrack(const reco::Conversion &c1, const reco::Conversion &c2)
const reco::PFCandidateCollection OniaPhotonConversionProducer::selectPFPhotons ( const reco::PFCandidateCollection pfcandidates)
private

Definition at line 320 of file OniaPhotonConversionProducer.cc.

References reco::PFCandidate::gamma.

320  {
321  reco::PFCandidateCollection pfphotons;
322  for (reco::PFCandidateCollection::const_iterator cand = pfcandidates.begin(); cand != pfcandidates.end(); ++cand){
323  if (cand->particleId() == reco::PFCandidate::gamma) pfphotons.push_back(*cand);
324  }
325  return pfphotons;
326 }
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates

Member Data Documentation

double OniaPhotonConversionProducer::_minDistanceOfApproachMaxCut
private

Definition at line 67 of file OniaPhotonConversionProducer.h.

double OniaPhotonConversionProducer::_minDistanceOfApproachMinCut
private

Definition at line 66 of file OniaPhotonConversionProducer.h.

double OniaPhotonConversionProducer::_trackchi2Cut
private

Definition at line 65 of file OniaPhotonConversionProducer.h.

double OniaPhotonConversionProducer::_vertexChi2ProbCut
private

Definition at line 64 of file OniaPhotonConversionProducer.h.

int OniaPhotonConversionProducer::convAlgo_
private

Definition at line 73 of file OniaPhotonConversionProducer.h.

edm::EDGetTokenT<reco::ConversionCollection> OniaPhotonConversionProducer::convCollectionToken_
private

Definition at line 55 of file OniaPhotonConversionProducer.h.

std::vector<int> OniaPhotonConversionProducer::convQuality_
private

Definition at line 74 of file OniaPhotonConversionProducer.h.

std::string OniaPhotonConversionProducer::convSelectionCuts_
private

Definition at line 76 of file OniaPhotonConversionProducer.h.

edm::EDGetTokenT<reco::PFCandidateCollection> OniaPhotonConversionProducer::pfCandidateCollectionToken_
private

Definition at line 57 of file OniaPhotonConversionProducer.h.

std::vector<double> OniaPhotonConversionProducer::pi0LargeWindow_
private

Definition at line 71 of file OniaPhotonConversionProducer.h.

bool OniaPhotonConversionProducer::pi0OnlineSwitch_
private

Definition at line 68 of file OniaPhotonConversionProducer.h.

std::vector<double> OniaPhotonConversionProducer::pi0SmallWindow_
private

Definition at line 70 of file OniaPhotonConversionProducer.h.

uint32_t OniaPhotonConversionProducer::sigmaTkVtxComp_
private

Definition at line 60 of file OniaPhotonConversionProducer.h.

edm::EDGetTokenT<reco::VertexCollection> OniaPhotonConversionProducer::thePVsToken_
private

Definition at line 56 of file OniaPhotonConversionProducer.h.

uint32_t OniaPhotonConversionProducer::TkMinNumOfDOF_
private

Definition at line 62 of file OniaPhotonConversionProducer.h.

bool OniaPhotonConversionProducer::wantCompatibleInnerHits_
private

Definition at line 61 of file OniaPhotonConversionProducer.h.

bool OniaPhotonConversionProducer::wantHighpurity_
private

Definition at line 63 of file OniaPhotonConversionProducer.h.

bool OniaPhotonConversionProducer::wantTkVtxCompatibility_
private

Definition at line 59 of file OniaPhotonConversionProducer.h.