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<>

Public Member Functions

 OniaPhotonConversionProducer (const edm::ParameterSet &ps)
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
bool hasAbilityToProduceInLumis () const final
 
bool hasAbilityToProduceInRuns () const final
 

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 &)
 
void endStream () override
 
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)
 
void produce (edm::Event &event, const edm::EventSetup &esetup) override
 
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::unique_ptr< StringCutObjectSelector< reco::Conversion > > convSelection_
 
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
 

Detailed Description

Select photon conversions and produce a conversion candidate collection

Definition at line 37 of file OniaPhotonConversionProducer.h.

Constructor & Destructor Documentation

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

Definition at line 53 of file OniaPhotonConversionProducer.cc.

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

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

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

325  {
326  // 2 windows are defined for Pi0 rejection, Conversions that, paired with others photons from the event, have an
327  // invariant mass inside the "small" window will be pizero_rejected and those that falls in the large window will
328  // be CheckPi0.
329  bool check_small = false;
330  bool check_large = false;
331 
332  float small1 = pi0SmallWindow_[0];
333  float small2 = pi0SmallWindow_[1];
334  float large1 = pi0LargeWindow_[0];
335  float large2 = pi0LargeWindow_[1];
336  for (reco::PFCandidateCollection::const_iterator photon = photons.begin(); photon!=photons.end(); ++photon) {
337  float inv = (conv.refittedPair4Momentum() + photon->p4()).M();
338  if (inv > large1 && inv < large2) {
339  check_large = true;
340  if (inv > small1 && inv < small2) {
341  check_small = true;
342  break;
343  }
344  }
345  }
346  pizero_rejected = check_small;
347  return check_large;
348 }
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 228 of file OniaPhotonConversionProducer.cc.

References begin, KineDebug3::count(), relativeConstraints::empty, end, training_settings::idx, lt_(), edm::second(), findQualityFiles::size, mathSSE::sqrt(), reco::Conversion::tracks(), and badGlobalMuonTaggersAOD_cff::vtx.

228  {
229  std::vector< std::pair< double, short> > idx[2];
230  short ik=-1;
231  for(auto const& tk : conv.tracks()){
232  ik++;
233  short count=-1;
234  for(auto const& vtx : priVtxs){
235  count++;
236 
237  double dz_= tk->dz(vtx.position());
238  double dzError_=tk->dzError();
239  dzError_=sqrt(dzError_*dzError_+vtx.covariance(2,2));
240  if(fabs(dz_)/dzError_ > sigmaTkVtxComp_) continue;
241  idx[ik].push_back(std::pair<double,short>(fabs(dz_),count));
242  }
243  if (idx[ik].empty()) return false;
244  if (idx[ik].size()>1) std::stable_sort(idx[ik].begin(),idx[ik].end(),lt_);
245  }
246  if (ik!=1) return false;
247  if (idx[0][0].second==idx[1][0].second) return true;
248  if (idx[0].size()>1 && idx[0][1].second==idx[1][0].second) return true;
249  if (idx[1].size()>1 && idx[0][0].second==idx[1][1].second) return true;
250 
251  return false;
252 }
size
Write out results.
U second(std::pair< T, U > const &p)
T sqrt(T t)
Definition: SSEVec.h:18
#define end
Definition: vmac.h:39
bool lt_(std::pair< double, short > a, std::pair< double, short > b)
#define begin
Definition: vmac.h:32
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 350 of file OniaPhotonConversionProducer.cc.

350  {
351  return reco::Candidate::LorentzVector(v.x(),v.y(), v.z(), v.t());
352 }
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
void OniaPhotonConversionProducer::endStream ( )
overrideprivate

Definition at line 355 of file OniaPhotonConversionProducer.cc.

References DEFINE_FWK_MODULE.

355  {
356 }
bool OniaPhotonConversionProducer::foundCompatibleInnerHits ( const reco::HitPattern hitPatA,
const reco::HitPattern hitPatB 
)
private

Definition at line 254 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::numberOfAllHits(), reco::HitPattern::trackerHitFilter(), and reco::HitPattern::validHitFilter().

254  {
255  size_t count=0;
256  uint32_t oldSubStr=0;
257  for (int i=0; i<hitPatA.numberOfAllHits(reco::HitPattern::HitCategory::TRACK_HITS) && count<2; i++) {
258  uint32_t hitA = hitPatA.getHitPattern(reco::HitPattern::HitCategory::TRACK_HITS,i);
259  if (!hitPatA.validHitFilter(hitA) || !hitPatA.trackerHitFilter(hitA)) continue;
260 
261  if(hitPatA.getSubStructure(hitA)==oldSubStr && hitPatA.getLayer(hitA)==oldSubStr)
262  continue;
263 
264  if(hitPatB.getTrackerMonoStereo(reco::HitPattern::HitCategory::TRACK_HITS,hitPatA.getSubStructure(hitA),hitPatA.getLayer(hitA)) != 0)
265  return true;
266 
267  oldSubStr=hitPatA.getSubStructure(hitA);
268  count++;
269  }
270  return false;
271 }
static uint32_t getLayer(uint16_t pattern)
Definition: HitPattern.h:759
uint16_t getTrackerMonoStereo(HitCategory category, uint16_t substr, uint16_t layer) const
Definition: HitPattern.cc:499
static bool validHitFilter(uint16_t pattern)
Definition: HitPattern.h:855
int numberOfAllHits(HitCategory category) const
Definition: HitPattern.h:875
static uint32_t getSubStructure(uint16_t pattern)
Definition: HitPattern.h:750
static bool trackerHitFilter(uint16_t pattern)
Definition: HitPattern.h:715
uint16_t getHitPattern(HitCategory category, int position) const
Definition: HitPattern.h:553
bool OniaPhotonConversionProducer::HighpuritySubset ( const reco::Conversion conv,
const reco::VertexCollection priVtxs 
)
private

Definition at line 274 of file OniaPhotonConversionProducer.cc.

References reco::Vertex::chi2(), ChiSquaredProbability(), reco::Conversion::conversionVertex(), reco::Conversion::distOfMinimumApproach(), mps_fire::i, reco::Vertex::ndof(), reco::Conversion::tracks(), badGlobalMuonTaggersAOD_cff::vtx, and reco::Conversion::zOfPrimaryVertexFromTracks().

Referenced by foundCompatibleInnerHits().

274  {
275  // select high purity conversions our way:
276  // vertex chi2 cut
278 
279  // d0 cut
280  // Find closest primary vertex
281  int closest_pv_index = 0;
282  int i=0;
283  for(auto const& vtx : priVtxs){
284  if( conv.zOfPrimaryVertexFromTracks( vtx.position() ) < conv.zOfPrimaryVertexFromTracks( priVtxs[closest_pv_index].position() ) ) closest_pv_index = i;
285  i++;
286  }
287  // Now check impact parameter wtr with the just found closest primary vertex
288  for(auto const& tk : conv.tracks()) if(-tk->dxy(priVtxs[closest_pv_index].position())*tk->charge()/tk->dxyError()<0) return false;
289 
290  // chi2 of single tracks
291  for(auto const& tk : conv.tracks()) if(tk->normalizedChi2() > _trackchi2Cut) return false;
292 
293  // dof for each track
294  for(auto const& tk : conv.tracks()) if(tk->ndof()< TkMinNumOfDOF_) return false;
295 
296  // distance of approach cut
298 
299  return true;
300 }
const reco::Vertex & conversionVertex() const
returns the reco conversion vertex
Definition: Conversion.h:98
double zOfPrimaryVertexFromTracks(const math::XYZPoint &myBeamSpot=math::XYZPoint()) const
Definition: Conversion.h:146
double distOfMinimumApproach() const
Definition: Conversion.h:126
double chi2() const
chi-squares
Definition: Vertex.h:98
float ChiSquaredProbability(double chiSquared, double nrDOF)
double ndof() const
Definition: Vertex.h:105
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 302 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().

302  {
303 
305  photonCand->setP4(convertVector(conv.refittedPair4Momentum()));
306  photonCand->setVertex(conv.conversionVertex().position());
307 
308  photonCand->addUserData<reco::Track>( "track0", *conv.tracks()[0] );
309  photonCand->addUserData<reco::Track>( "track1", *conv.tracks()[1] );
310 
311  return photonCand;
312 
313 }
const reco::Vertex & conversionVertex() const
returns the reco conversion vertex
Definition: Conversion.h:98
Analysis-level particle class.
reco::Candidate::LorentzVector convertVector(const math::XYZTLorentzVectorF &)
const Point & position() const
position
Definition: Vertex.h:109
void setVertex(const Point &vertex) override
set vertex
math::XYZTLorentzVectorF refittedPair4Momentum() const
Conversion track pair 4-momentum from the tracks refitted with vertex constraint. ...
Definition: Conversion.cc:235
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 177 of file OniaPhotonConversionProducer.cc.

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

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

Definition at line 85 of file OniaPhotonConversionProducer.cc.

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

85  {
86 
87  std::unique_ptr<reco::ConversionCollection> outCollection(new reco::ConversionCollection);
88  std::unique_ptr<pat::CompositeCandidateCollection> patoutCollection(new pat::CompositeCandidateCollection);
89 
91  event.getByToken(thePVsToken_, priVtxs);
92 
94  event.getByToken(convCollectionToken_,pConv);
95 
97  event.getByToken(pfCandidateCollectionToken_,pfcandidates);
98 
99  const reco::PFCandidateCollection pfphotons = selectPFPhotons(*pfcandidates);
100 
101  for(reco::ConversionCollection::const_iterator conv = pConv->begin(); conv != pConv->end(); ++conv){
102 
103  if (! ( *convSelection_)(*conv)){
104  continue; // selection string
105  }
106  if (convAlgo_ != 0 && conv->algo()!= convAlgo_){
107  continue; // select algorithm
108  }
109  if(!convQuality_.empty()){
110  bool flagsok=true;
111  for (std::vector<int>::const_iterator flag = convQuality_.begin(); flag!=convQuality_.end(); ++flag){
113  if (!conv->quality(q)) {
114  flagsok=false;
115  break;
116  }
117  }
118  if (!flagsok){
119  continue;
120  }
121  }
122  outCollection->push_back(*conv);
123  }
124 
125  removeDuplicates(*outCollection);
126 
127  for (reco::ConversionCollection::const_iterator conv = outCollection->begin(); conv != outCollection->end(); ++conv){
128 
129  bool flag1 = true;
130  bool flag2 = true;
131  bool flag3 = true;
132  bool flag4 = true;
133 
134  // The logic implies that by default this flags are true and if the check are not wanted conversions are saved.
135  // If checks are required and failed then don't save the conversion.
136 
137  bool flagTkVtxCompatibility = true;
138  if (! checkTkVtxCompatibility(*conv,*priVtxs.product())) {
139  flagTkVtxCompatibility = false;
141  flag1 = false;
142  }
143  }
144  bool flagCompatibleInnerHits = false;
145  if (conv->tracks().size()==2) {
146  reco::HitPattern hitPatA=conv->tracks().at(0)->hitPattern();
147  reco::HitPattern hitPatB=conv->tracks().at(1)->hitPattern();
148  if ( foundCompatibleInnerHits(hitPatA,hitPatB) && foundCompatibleInnerHits(hitPatB,hitPatA) ) flagCompatibleInnerHits = true;
149  }
150  if (wantCompatibleInnerHits_ && ! flagCompatibleInnerHits) {
151  flag2 = false;
152  }
153  bool flagHighpurity = true;
154  if (!HighpuritySubset(*conv,*priVtxs.product())) {
155  flagHighpurity = false;
156  if (wantHighpurity_) {
157  flag3 = false;
158  }
159  }
160  bool pizero_rejected = false;
161  bool large_pizero_window = CheckPi0(*conv, pfphotons, pizero_rejected);
162  if (pi0OnlineSwitch_ && pizero_rejected) {
163  flag4 = false;
164  }
165 
166  int flags = 0;
167  if (flag1 && flag2 && flag3 && flag4){
168  flags = PackFlags(*conv,flagTkVtxCompatibility,flagCompatibleInnerHits,flagHighpurity,pizero_rejected,large_pizero_window);
169  std::unique_ptr<pat::CompositeCandidate> pat_conv(makePhotonCandidate(*conv));
170  pat_conv->addUserInt("flags",flags);
171  patoutCollection->push_back(*pat_conv);
172  }
173  }
174  event.put(std::move(patoutCollection),"conversions");
175 }
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
std::unique_ptr< StringCutObjectSelector< reco::Conversion > > convSelection_
bool checkTkVtxCompatibility(const reco::Conversion &, const reco::VertexCollection &)
int PackFlags(const reco::Conversion &, bool, bool, bool, bool, bool)
edm::EDGetTokenT< reco::VertexCollection > thePVsToken_
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:74
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:511
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 210 of file OniaPhotonConversionProducer.cc.

References ConversionEqualByTrack(), and ConversionLessByChi2().

210  {
211  // first sort from high to low chi2 prob
212  std::sort(c.begin(),c.end(),ConversionLessByChi2);
213  int iter1 = 0;
214  // Cycle over all the elements of the collection and compare to all the following,
215  // if two elements have at least one track in common delete the element with the lower chi2
216  while(iter1 < (( (int) c.size() ) - 1) ){
217  int iter2 = iter1+1;
218  while( iter2 < (int) c.size()) if(ConversionEqualByTrack( c[iter1], c[iter2] ) ){
219  c.erase( c.begin() + iter2 );
220  }else{
221  iter2++; // Increment index only if this element is no duplicate.
222  // If it is duplicate check again the same index since the vector rearranged elements index after erasing
223  }
224  iter1++;
225  }
226 }
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 316 of file OniaPhotonConversionProducer.cc.

References reco::PFCandidate::gamma.

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

Member Data Documentation

double OniaPhotonConversionProducer::_minDistanceOfApproachMaxCut
private

Definition at line 68 of file OniaPhotonConversionProducer.h.

double OniaPhotonConversionProducer::_minDistanceOfApproachMinCut
private

Definition at line 67 of file OniaPhotonConversionProducer.h.

double OniaPhotonConversionProducer::_trackchi2Cut
private

Definition at line 66 of file OniaPhotonConversionProducer.h.

double OniaPhotonConversionProducer::_vertexChi2ProbCut
private

Definition at line 65 of file OniaPhotonConversionProducer.h.

int OniaPhotonConversionProducer::convAlgo_
private

Definition at line 74 of file OniaPhotonConversionProducer.h.

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

Definition at line 56 of file OniaPhotonConversionProducer.h.

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

Definition at line 75 of file OniaPhotonConversionProducer.h.

std::unique_ptr<StringCutObjectSelector<reco::Conversion> > OniaPhotonConversionProducer::convSelection_
private

Definition at line 78 of file OniaPhotonConversionProducer.h.

std::string OniaPhotonConversionProducer::convSelectionCuts_
private

Definition at line 77 of file OniaPhotonConversionProducer.h.

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

Definition at line 58 of file OniaPhotonConversionProducer.h.

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

Definition at line 72 of file OniaPhotonConversionProducer.h.

bool OniaPhotonConversionProducer::pi0OnlineSwitch_
private

Definition at line 69 of file OniaPhotonConversionProducer.h.

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

Definition at line 71 of file OniaPhotonConversionProducer.h.

uint32_t OniaPhotonConversionProducer::sigmaTkVtxComp_
private

Definition at line 61 of file OniaPhotonConversionProducer.h.

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

Definition at line 57 of file OniaPhotonConversionProducer.h.

uint32_t OniaPhotonConversionProducer::TkMinNumOfDOF_
private

Definition at line 63 of file OniaPhotonConversionProducer.h.

bool OniaPhotonConversionProducer::wantCompatibleInnerHits_
private

Definition at line 62 of file OniaPhotonConversionProducer.h.

bool OniaPhotonConversionProducer::wantHighpurity_
private

Definition at line 64 of file OniaPhotonConversionProducer.h.

bool OniaPhotonConversionProducer::wantTkVtxCompatibility_
private

Definition at line 60 of file OniaPhotonConversionProducer.h.