CMS 3D CMS Logo

OniaPhotonConversionProducer.cc
Go to the documentation of this file.
8 
12 
15 
17 
18 #include <TMath.h>
19 #include <boost/foreach.hpp>
20 #include <vector>
21 
22 
23 // to order from high to low ProbChi2
25  return TMath::Prob(c1.conversionVertex().chi2(),c1.conversionVertex().ndof()) > TMath::Prob(c2.conversionVertex().chi2(),c2.conversionVertex().ndof());
26 }
27 
28 
30  bool atLeastOneInCommon=false;
31  BOOST_FOREACH(const edm::RefToBase<reco::Track> tk1, c1.tracks()){
32  BOOST_FOREACH(const edm::RefToBase<reco::Track> tk2, c2.tracks()){
33  if(tk1==tk2){
34  atLeastOneInCommon=true;
35  break;
36  }
37  }
38  }
39  return atLeastOneInCommon;
40 }
41 
42 bool lt_(std::pair<double,short> a, std::pair<double,short> b) {
43  return a.first < b.first; }
44 
45 // define operator== for conversions, those with at least one track in common
46 namespace reco {
47  bool operator==(const reco::Conversion& c1, const reco::Conversion& c2) {
48  return c1.tracks()[0] == c2.tracks()[0] ||
49  c1.tracks()[1] == c2.tracks()[1] ||
50  c1.tracks()[1] == c2.tracks()[0] ||
51  c1.tracks()[0] == c2.tracks()[1] ;
52  }
53 }
54 
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 }
85 
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;
143  if (wantTkVtxCompatibility_) {
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);
172  std::unique_ptr<pat::CompositeCandidate> pat_conv(makePhotonCandidate(*conv));
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 }
181 
182 int OniaPhotonConversionProducer::PackFlags(const reco::Conversion& conv, bool flagTkVtxCompatibility,
183  bool flagCompatibleInnerHits, bool flagHighpurity,
184  bool pizero_rejected, bool large_pizero_window ) {
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 }
211 
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 }
232 
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  if(fabs(dz_)/dzError_ > sigmaTkVtxComp_) continue;
246  idx[ik].push_back(std::pair<double,short>(fabs(dz_),count));
247  }
248  if (idx[ik].size()==0) return false;
249  if (idx[ik].size()>1) std::stable_sort(idx[ik].begin(),idx[ik].end(),lt_);
250  }
251  if (ik!=1) return false;
252  if (idx[0][0].second==idx[1][0].second) return true;
253  if (idx[0].size()>1 && idx[0][1].second==idx[1][0].second) return true;
254  if (idx[1].size()>1 && idx[0][0].second==idx[1][1].second) return true;
255 
256  return false;
257 }
258 
260  size_t count=0;
261  uint32_t oldSubStr=0;
262  for (int i=0; i<hitPatA.numberOfHits(reco::HitPattern::HitCategory::TRACK_HITS) && count<2; i++) {
263  uint32_t hitA = hitPatA.getHitPattern(reco::HitPattern::HitCategory::TRACK_HITS,i);
264  if (!hitPatA.validHitFilter(hitA) || !hitPatA.trackerHitFilter(hitA)) continue;
265 
266  if(hitPatA.getSubStructure(hitA)==oldSubStr && hitPatA.getLayer(hitA)==oldSubStr)
267  continue;
268 
269  if(hitPatB.getTrackerMonoStereo(reco::HitPattern::HitCategory::TRACK_HITS,hitPatA.getSubStructure(hitA),hitPatA.getLayer(hitA)) != 0)
270  return true;
271 
272  oldSubStr=hitPatA.getSubStructure(hitA);
273  count++;
274  }
275  return false;
276 }
277 
280  // select high purity conversions our way:
281  // vertex chi2 cut
282  if(ChiSquaredProbability(conv.conversionVertex().chi2(),conv.conversionVertex().ndof())< _vertexChi2ProbCut) return false;
283 
284  // d0 cut
285  // Find closest primary vertex
286  int closest_pv_index = 0;
287  int i=0;
288  BOOST_FOREACH(const reco::Vertex& vtx,priVtxs){
289  if( conv.zOfPrimaryVertexFromTracks( vtx.position() ) < conv.zOfPrimaryVertexFromTracks( priVtxs[closest_pv_index].position() ) ) closest_pv_index = i;
290  i++;
291  }
292  // Now check impact parameter wtr with the just found closest primary vertex
293  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;
294 
295  // chi2 of single tracks
296  BOOST_FOREACH(const edm::RefToBase<reco::Track> tk, conv.tracks()) if(tk->normalizedChi2() > _trackchi2Cut) return false;
297 
298  // dof for each track
299  BOOST_FOREACH(const edm::RefToBase<reco::Track> tk, conv.tracks()) if(tk->ndof()< TkMinNumOfDOF_) return false;
300 
301  // distance of approach cut
302  if (conv.distOfMinimumApproach() < _minDistanceOfApproachMinCut || conv.distOfMinimumApproach() > _minDistanceOfApproachMaxCut ) return false;
303 
304  return true;
305 }
306 
308 
310  photonCand->setP4(convertVector(conv.refittedPair4Momentum()));
311  photonCand->setVertex(conv.conversionVertex().position());
312 
313  photonCand->addUserData<reco::Track>( "track0", *conv.tracks()[0] );
314  photonCand->addUserData<reco::Track>( "track1", *conv.tracks()[1] );
315 
316  return photonCand;
317 
318 }
319 
320 // create a collection of PF photons
322  reco::PFCandidateCollection pfphotons;
323  for (reco::PFCandidateCollection::const_iterator cand = pfcandidates.begin(); cand != pfcandidates.end(); ++cand){
324  if (cand->particleId() == reco::PFCandidate::gamma) pfphotons.push_back(*cand);
325  }
326  return pfphotons;
327 }
328 
330  bool &pizero_rejected ) {
331  // 2 windows are defined for Pi0 rejection, Conversions that, paired with others photons from the event, have an
332  // invariant mass inside the "small" window will be pizero_rejected and those that falls in the large window will
333  // be CheckPi0.
334  bool check_small = false;
335  bool check_large = false;
336 
337  float small1 = pi0SmallWindow_[0];
338  float small2 = pi0SmallWindow_[1];
339  float large1 = pi0LargeWindow_[0];
340  float large2 = pi0LargeWindow_[1];
341  for (reco::PFCandidateCollection::const_iterator photon = photons.begin(); photon!=photons.end(); ++photon) {
342  float inv = (conv.refittedPair4Momentum() + photon->p4()).M();
343  if (inv > large1 && inv < large2) {
344  check_large = true;
345  if (inv > small1 && inv < small2) {
346  check_small = true;
347  break;
348  }
349  }
350  }
351  pizero_rejected = check_small;
352  return check_large;
353 }
354 
356  return reco::Candidate::LorentzVector(v.x(),v.y(), v.z(), v.t());
357 }
358 
359 
361 }
362 //define this as a plug-in
size
Write out results.
const reco::Vertex & conversionVertex() const
returns the reco conversion vertex
Definition: Conversion.h:97
Analysis-level particle class.
T getParameter(std::string const &) const
static uint32_t getLayer(uint16_t pattern)
Definition: HitPattern.h:700
bool foundCompatibleInnerHits(const reco::HitPattern &hitPatA, const reco::HitPattern &hitPatB)
ConversionAlgorithm algo() const
Definition: Conversion.h:223
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:561
static HepMC::IO_HEPEVT conv
double dxyError() const
error on dxy
Definition: TrackBase.h:796
bool quality(ConversionQuality q) const
Definition: Conversion.h:181
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
reco::Candidate::LorentzVector convertVector(const math::XYZTLorentzVectorF &)
uint16_t getTrackerMonoStereo(HitCategory category, uint16_t substr, uint16_t layer) const
Definition: HitPattern.cc:454
double zOfPrimaryVertexFromTracks(const math::XYZPoint &myBeamSpot=math::XYZPoint()) const
Definition: Conversion.h:145
double distOfMinimumApproach() const
Definition: Conversion.h:125
bool HighpuritySubset(const reco::Conversion &, const reco::VertexCollection &)
std::vector< Variable::Flags > flags
Definition: MVATrainer.cc:135
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
bool ConversionLessByChi2(const reco::Conversion &c1, const reco::Conversion &c2)
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
std::vector< Conversion > ConversionCollection
collectin of Conversion objects
Definition: ConversionFwd.h:9
static bool validHitFilter(uint16_t pattern)
Definition: HitPattern.h:787
U second(std::pair< T, U > const &p)
virtual void produce(edm::Event &event, const edm::EventSetup &esetup)
bool operator==(const reco::Conversion &c1, const reco::Conversion &c2)
bool ConversionEqualByTrack(const reco::Conversion &c1, const reco::Conversion &c2)
bool checkTkVtxCompatibility(const reco::Conversion &, const reco::VertexCollection &)
double ndof() const
number of degrees of freedom of the fit
Definition: TrackBase.h:555
T sqrt(T t)
Definition: SSEVec.h:18
int PackFlags(const reco::Conversion &, bool, bool, bool, bool, bool)
double chi2() const
chi-squares
Definition: Vertex.h:98
float ChiSquaredProbability(double chiSquared, double nrDOF)
#define end
Definition: vmac.h:37
virtual void setVertex(const Point &vertex)
set vertex
static uint32_t getSubStructure(uint16_t pattern)
Definition: HitPattern.h:691
double ndof() const
Definition: Vertex.h:105
const reco::PFCandidateCollection selectPFPhotons(const reco::PFCandidateCollection &)
OniaPhotonConversionProducer(const edm::ParameterSet &ps)
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:609
double dzError() const
error on dz
Definition: TrackBase.h:814
math::XYZTLorentzVectorF refittedPair4Momentum() const
Conversion track pair 4-momentum from the tracks refitted with vertex constraint. ...
Definition: Conversion.cc:235
bool lt_(std::pair< double, short > a, std::pair< double, short > b)
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
static bool trackerHitFilter(uint16_t pattern)
Definition: HitPattern.h:677
std::vector< CompositeCandidate > CompositeCandidateCollection
T const * product() const
Definition: Handle.h:81
double b
Definition: hdecay.h:120
virtual void setP4(const LorentzVector &p4) final
set 4-momentum
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
fixed size matrix
#define begin
Definition: vmac.h:30
bool CheckPi0(const reco::Conversion &, const reco::PFCandidateCollection &, bool &)
pat::CompositeCandidate * makePhotonCandidate(const reco::Conversion &)
double a
Definition: hdecay.h:121
int charge() const
track electric charge
Definition: TrackBase.h:567
void removeDuplicates(reco::ConversionCollection &)
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< float > > XYZTLorentzVectorF
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:22
uint16_t getHitPattern(HitCategory category, int position) const
Definition: HitPattern.h:515
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:591
void addUserData(const std::string &label, const T &data, bool transientOnly=false, bool overwrite=false)
Definition: PATObject.h:309
def move(src, dest)
Definition: eostools.py:510
std::vector< edm::RefToBase< reco::Track > > const & tracks() const
vector of track to base references
Definition: Conversion.cc:176
Definition: event.py:1
int numberOfHits(HitCategory category) const
Definition: HitPattern.h:807