CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Classes | Public Member Functions | Private Member Functions | Private Attributes | Friends
BtoCharmDecayVertexMerger Class Reference
Inheritance diagram for BtoCharmDecayVertexMerger:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

Classes

struct  vertexProxy
 

Public Member Functions

 BtoCharmDecayVertexMerger (const edm::ParameterSet &params)
 
virtual void produce (edm::Event &event, const edm::EventSetup &es)
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
boost::function< void(const
BranchDescription &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 

Private Member Functions

GlobalVector flightDirection (reco::Vertex &pv, reco::Vertex &sv)
 
void resolveBtoDchain (std::vector< vertexProxy > &coll, unsigned int i, unsigned int k)
 

Private Attributes

double maxPtreltomerge
 
double minCosPAtomerge
 
double minDRForUnique
 
edm::InputTag primaryVertexCollection
 
reco::Vertex pv
 
edm::InputTag secondaryVertexCollection
 
double vecSumIMCUTForUnique
 

Friends

bool operator< (vertexProxy v1, vertexProxy v2)
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
typedef WorkerT< EDProducerWorkerType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDProducer
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::ProducerBase
template<class TProducer , class TMethod >
void callWhenNewProductsRegistered (TProducer *iProd, TMethod iMethod)
 

Detailed Description

Definition at line 23 of file BtoCharmDecayVertexMerger.cc.

Constructor & Destructor Documentation

BtoCharmDecayVertexMerger::BtoCharmDecayVertexMerger ( const edm::ParameterSet params)

Definition at line 60 of file BtoCharmDecayVertexMerger.cc.

60  :
61  primaryVertexCollection(params.getParameter<edm::InputTag>("primaryVertices")),
62  secondaryVertexCollection(params.getParameter<edm::InputTag>("secondaryVertices")),
63  minDRForUnique(params.getUntrackedParameter<double>("minDRUnique",0.4)),
64  vecSumIMCUTForUnique(params.getUntrackedParameter<double>("minvecSumIMifsmallDRUnique",5.5)),
65  minCosPAtomerge(params.getUntrackedParameter<double>("minCosPAtomerge",0.99)),
66  maxPtreltomerge(params.getUntrackedParameter<double>("maxPtreltomerge",7777.0))
67  // maxFraction(params.getParameter<double>("maxFraction")),
68  // minSignificance(params.getParameter<double>("minSignificance"))
69 {
70  produces<reco::VertexCollection>();
71 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const

Member Function Documentation

GlobalVector BtoCharmDecayVertexMerger::flightDirection ( reco::Vertex pv,
reco::Vertex sv 
)
private

Definition at line 235 of file BtoCharmDecayVertexMerger.cc.

References reco::Vertex::position().

Referenced by resolveBtoDchain().

235  {
236  GlobalVector res(sv.position().X() - pv.position().X(),
237  sv.position().Y() - pv.position().Y(),
238  sv.position().Z() - pv.position().Z());
239  return res;
240 }
const Point & position() const
position
Definition: Vertex.h:93
void BtoCharmDecayVertexMerger::produce ( edm::Event event,
const edm::EventSetup es 
)
virtual

Implements edm::EDProducer.

Definition at line 73 of file BtoCharmDecayVertexMerger.cc.

References edm::Event::getByLabel(), primaryVertexCollection, edm::Event::put(), pv, dt_dqm_sourceclient_common_cff::reco, resolveBtoDchain(), secondaryVertexCollection, and python.multivaluedict::sort().

73  {
74 
75  using namespace reco;
76  // PV
78  iEvent.getByLabel(primaryVertexCollection, PVcoll);
79 
80  if(PVcoll->size()!=0) {
81 
82  const reco::VertexCollection pvc = *( PVcoll.product());
83  pv = pvc[0];
84 
85  // get the IVF collection
86  edm::Handle<VertexCollection> secondaryVertices;
87  iEvent.getByLabel(secondaryVertexCollection, secondaryVertices);
88 
89 
90 
91  //loop over vertices, fill into own collection for sorting
92  std::vector<vertexProxy> vertexProxyColl;
93  for(std::vector<reco::Vertex>::const_iterator sv = secondaryVertices->begin();
94  sv != secondaryVertices->end(); ++sv) {
95  vertexProxy vtx = {*sv,(*sv).p4().M(),(*sv).tracksSize()};
96  vertexProxyColl.push_back( vtx );
97  }
98 
99  // sort the vertices by mass and track multiplicity
100  sort( vertexProxyColl.begin(), vertexProxyColl.end());
101 
102 
103  // loop forward over all vertices
104  for(unsigned int iVtx=0; iVtx < vertexProxyColl.size(); iVtx++){
105 
106  // nested loop backwards (in order to start from light masses)
107  // check all vertices against each other for B->D chain
108  for(unsigned int kVtx=vertexProxyColl.size()-1; kVtx>iVtx; kVtx--){
109  // remove D vertices from the collection and add the tracks to the original one
110  resolveBtoDchain(vertexProxyColl, iVtx, kVtx);
111  }
112  }
113 
114  // now create new vertex collection and add to event
115  VertexCollection *bvertices = new VertexCollection();
116  for(std::vector<vertexProxy>::iterator it=vertexProxyColl.begin(); it!=vertexProxyColl.end(); it++) bvertices->push_back((*it).vert);
117  std::auto_ptr<VertexCollection> bvertColl(bvertices);
118  iEvent.put(bvertColl);
119  }
120  else{
121  std::auto_ptr<VertexCollection> bvertCollEmpty(new VertexCollection);
122  iEvent.put(bvertCollEmpty);
123  }
124 }
void resolveBtoDchain(std::vector< vertexProxy > &coll, unsigned int i, unsigned int k)
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
int iEvent
Definition: GenABIO.cc:243
void BtoCharmDecayVertexMerger::resolveBtoDchain ( std::vector< vertexProxy > &  coll,
unsigned int  i,
unsigned int  k 
)
private

Definition at line 128 of file BtoCharmDecayVertexMerger.cc.

References Geom::deltaR(), reco::SecondaryVertex::dist3d(), dot(), spr::find(), flightDirection(), i, gen::k, maxPtreltomerge, minCosPAtomerge, minDRForUnique, edm::RefVector< C, T, F >::push_back(), pv, dt_dqm_sourceclient_common_cff::reco, reco::Vertex::refittedTrack(), mathSSE::sqrt(), lumiQTWidget::t, groupFilesInBlocks::temp, tracks_, reco::Vertex::tracks_begin(), reco::Vertex::tracks_end(), reco::Vertex::trackWeight(), Measurement1D::value(), vecSumIMCUTForUnique, and CommonMethods::weight().

Referenced by produce().

128  {
129  using namespace reco;
130 
131  GlobalVector momentum1 = GlobalVector(coll[i].vert.p4().X(), coll[i].vert.p4().Y(), coll[i].vert.p4().Z());
132  GlobalVector momentum2 = GlobalVector(coll[k].vert.p4().X(), coll[k].vert.p4().Y(), coll[k].vert.p4().Z());
133 
134  reco::SecondaryVertex sv1(pv, coll[i].vert, momentum1 , true);
135  reco::SecondaryVertex sv2(pv, coll[k].vert, momentum2 , true);
136 
137 
138  // find out which one is near and far
139  reco::SecondaryVertex svNear = sv1;
140  reco::SecondaryVertex svFar = sv2;
141  GlobalVector momentumNear = momentum1;
142  GlobalVector momentumFar = momentum2;
143 
144  // swap if it is the other way around
145  if(sv1.dist3d().value() >= sv2.dist3d().value()){
146  svNear = sv2;
147  svFar = sv1;
148  momentumNear = momentum2;
149  momentumFar = momentum1;
150  }
151  GlobalVector nearToFar = flightDirection( svNear, svFar);
152  GlobalVector pvToNear = flightDirection( pv, svNear);
153 
154  double cosPA = nearToFar.dot(momentumFar) / momentumFar.mag()/ nearToFar.mag();
155  double cosa = pvToNear. dot(momentumFar) / pvToNear.mag() / momentumFar.mag();
156  double ptrel = sqrt(1.0 - cosa*cosa)* momentumFar.mag();
157 
158  double vertexDeltaR = Geom::deltaR(flightDirection(pv, sv1), flightDirection(pv, sv2) );
159 
160  // create a set of all tracks from both vertices, avoid double counting by using a std::set<>
161  std::set<reco::TrackRef> trackrefs;
162  // first vertex
163  for(reco::Vertex::trackRef_iterator ti = sv1.tracks_begin(); ti!=sv1.tracks_end(); ti++){
164  if(sv1.trackWeight(*ti)>0.5){
165  reco::TrackRef t = ti->castTo<reco::TrackRef>();
166  trackrefs.insert(t);
167  }
168  }
169  // second vertex
170  for(reco::Vertex::trackRef_iterator ti = sv2.tracks_begin(); ti!=sv2.tracks_end(); ti++){
171  if(sv2.trackWeight(*ti)>0.5){
172  reco::TrackRef t = ti->castTo<reco::TrackRef>();
173  trackrefs.insert(t);
174  }
175  }
176 
177  // now calculate one LorentzVector from the track momenta
178  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > mother;
179  for(std::set<reco::TrackRef>::const_iterator it = trackrefs.begin(); it!= trackrefs.end(); it++){
180  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > temp ( (*it)->px(),(*it)->py(),(*it)->pz(), 0.13957 );
181  mother += temp;
182  }
183 
184 
185  // check if criteria for merging are fulfilled
186  if(vertexDeltaR<minDRForUnique && mother.M()<vecSumIMCUTForUnique && cosPA>minCosPAtomerge && fabs(ptrel)<maxPtreltomerge ) {
187 
188 
189  // add tracks of second vertex which are missing to first vertex
190  // loop over the second
191  bool bFoundDuplicate=false;
192  for(reco::Vertex::trackRef_iterator ti = sv2.tracks_begin(); ti!=sv2.tracks_end(); ti++){
193  reco::Vertex::trackRef_iterator it = find(sv1.tracks_begin(), sv1.tracks_end(), *ti);
194  if (it==sv1.tracks_end()) coll[i].vert.add( *ti, sv2.refittedTrack(*ti), sv2.trackWeight(*ti) );
195  else bFoundDuplicate=true;
196  }
197  // in case a duplicate track is found, need to create the full track list in the vertex from scratch because we need to modify the weight.
198  // the weight must be the larger one, otherwise we may have outliers which are not real outliers
199 
200  if(bFoundDuplicate){
201  // create backup track containers from main vertex
202  std::vector<TrackBaseRef > tracks_;
203  std::vector<Track> refittedTracks_;
204  std::vector<float> weights_;
205  for(reco::Vertex::trackRef_iterator it = coll[i].vert.tracks_begin(); it!=coll[i].vert.tracks_end(); it++) {
206  tracks_.push_back( *it);
207  refittedTracks_.push_back( coll[i].vert.refittedTrack(*it));
208  weights_.push_back( coll[i].vert.trackWeight(*it) );
209  }
210  // delete tracks and add all tracks back, and check in which vertex the weight is larger
211  coll[i].vert.removeTracks();
212  std::vector<Track>::iterator it2 = refittedTracks_.begin();
213  std::vector<float>::iterator it3 = weights_.begin();
214  for(reco::Vertex::trackRef_iterator it = tracks_.begin(); it!=tracks_.end(); it++, it2++, it3++){
215  float weight = *it3;
216  float weight2= sv2.trackWeight(*it);
217  Track refittedTrackWithLargerWeight = *it2;
218  if( weight2 >weight) {
219  weight = weight2;
220  refittedTrackWithLargerWeight = sv2.refittedTrack(*it);
221  }
222  coll[i].vert.add(*it , refittedTrackWithLargerWeight , weight);
223  }
224  }
225 
226  // remove the second vertex from the collection
227  coll.erase( coll.begin() + k );
228  }
229 
230 }
int i
Definition: DBlmapReader.cc:9
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
GlobalVector flightDirection(reco::Vertex &pv, reco::Vertex &sv)
T sqrt(T t)
Definition: SSEVec.h:46
int k[5][pyjets_maxn]
JetCorrectorParametersCollection coll
Definition: classes.h:14
double deltaR(const Vector1 &v1, const Vector2 &v2)
Definition: VectorUtil.h:84
T dot(const Basic3DVector &v) const
Scalar product, or &quot;dot&quot; product, with a vector of same type.
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:38
void push_back(value_type const &ref)
Add a Ref&lt;C, T&gt; to the RefVector.
Definition: RefVector.h:64
const reco::PFCandidateRefVector & tracks_
Global3DVector GlobalVector
Definition: GlobalVector.h:10

Friends And Related Function Documentation

bool operator< ( vertexProxy  v1,
vertexProxy  v2 
)
friend

Definition at line 48 of file BtoCharmDecayVertexMerger.cc.

48  {
49  if(v1.ntracks>2 && v2.ntracks<3) return true;
50  if(v1.ntracks<3 && v2.ntracks>2) return false;
51  return (v1.invm>v2.invm);
52  }

Member Data Documentation

double BtoCharmDecayVertexMerger::maxPtreltomerge
private

Definition at line 38 of file BtoCharmDecayVertexMerger.cc.

Referenced by resolveBtoDchain().

double BtoCharmDecayVertexMerger::minCosPAtomerge
private

Definition at line 38 of file BtoCharmDecayVertexMerger.cc.

Referenced by resolveBtoDchain().

double BtoCharmDecayVertexMerger::minDRForUnique
private

Definition at line 38 of file BtoCharmDecayVertexMerger.cc.

Referenced by resolveBtoDchain().

edm::InputTag BtoCharmDecayVertexMerger::primaryVertexCollection
private

Definition at line 32 of file BtoCharmDecayVertexMerger.cc.

Referenced by produce().

reco::Vertex BtoCharmDecayVertexMerger::pv
private

Definition at line 34 of file BtoCharmDecayVertexMerger.cc.

Referenced by produce(), and resolveBtoDchain().

edm::InputTag BtoCharmDecayVertexMerger::secondaryVertexCollection
private

Definition at line 33 of file BtoCharmDecayVertexMerger.cc.

Referenced by produce().

double BtoCharmDecayVertexMerger::vecSumIMCUTForUnique
private

Definition at line 38 of file BtoCharmDecayVertexMerger.cc.

Referenced by resolveBtoDchain().