CMS 3D CMS Logo

BPHKinematicFit.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * \author Paolo Ronchese INFN Padova
5  *
6  */
7 
8 //-----------------------
9 // This Class' Header --
10 //-----------------------
12 
13 //-------------------------------
14 // Collaborating Class Headers --
15 //-------------------------------
26 
27 //---------------
28 // C++ Headers --
29 //---------------
30 #include <iostream>
31 
32 using namespace std;
33 
34 
35 //-------------------
36 // Initializations --
37 //-------------------
38 
39 
40 //----------------
41 // Constructors --
42 //----------------
44  BPHDecayVertex( 0 ),
45  massConst( -1.0 ),
46  massSigma( -1.0 ),
47  oldKPs( true ),
48  oldFit( true ),
49  oldMom( true ),
50  kinTree( 0 ) {
51 }
52 
53 
55  BPHDecayVertex( ptr, 0 ),
56  massConst( -1.0 ),
57  massSigma( -1.0 ),
58  oldKPs( true ),
59  oldFit( true ),
60  oldMom( true ),
61  kinTree( 0 ) {
62  map<const reco::Candidate*,const reco::Candidate*> iMap;
63  const vector<const reco::Candidate*>& daug = daughters();
64  const vector<Component>& list = ptr->componentList();
65  int i;
66  int n = daug.size();
67  for ( i = 0; i < n; ++i ) {
68  const reco::Candidate* cand = daug[i];
69  iMap[originalReco( cand )] = cand;
70  }
71  for ( i = 0; i < n; ++i ) {
72  const Component& c = list[i];
73  dMSig[iMap[c.cand]] = c.msig;
74  }
75  const vector<BPHRecoConstCandPtr>& dComp = daughComp();
76  int j;
77  int m = dComp.size();
78  for ( j = 0; j < m; ++j ) {
79  const map<const reco::Candidate*,double>& dMap = dComp[j]->dMSig;
80  dMSig.insert( dMap.begin(), dMap.end() );
81  }
82 }
83 
84 //--------------
85 // Destructor --
86 //--------------
88 }
89 
90 //--------------
91 // Operations --
92 //--------------
93 void BPHKinematicFit::setConstraint( double mass, double sigma ) {
94  oldFit = oldMom = true;
95  massConst = mass;
96  massSigma = sigma;
97  return;
98 }
99 
100 
102  return massConst;
103 }
104 
105 
107  return massSigma;
108 }
109 
110 
111 const vector<RefCountedKinematicParticle>& BPHKinematicFit::kinParticles()
112  const {
113  if ( oldKPs ) buildParticles();
114  return allParticles;
115 }
116 
117 
118 vector<RefCountedKinematicParticle> BPHKinematicFit::kinParticles(
119  const vector<string>& names ) const {
120  if ( oldKPs ) buildParticles();
121  const vector<const reco::Candidate*>& daugs = daughFull();
122  vector<RefCountedKinematicParticle> plist;
123  if ( allParticles.size() != daugs.size() ) return plist;
124  set<RefCountedKinematicParticle> pset;
125  int i;
126  int n = names.size();
127  int m = daugs.size();
128  plist.reserve( m );
129  for ( i = 0; i < n; ++i ) {
130  const string& pname = names[i];
131  if ( pname == "*" ) {
132  int j = m;
133  while ( j-- ) {
135  if ( pset.find( kp ) != pset.end() ) continue;
136  plist.push_back( kp );
137  pset .insert ( kp );
138  }
139  break;
140  }
141  map<const reco::Candidate*,
142  RefCountedKinematicParticle>::const_iterator iter = kinMap.find(
143  getDaug( pname ) );
144  map<const reco::Candidate*,
145  RefCountedKinematicParticle>::const_iterator iend = kinMap.end();
146  if ( iter != iend ) {
147  const RefCountedKinematicParticle& kp = iter->second;
148  if ( pset.find( kp ) != pset.end() ) continue;
149  plist.push_back( kp );
150  pset .insert ( kp );
151  }
152  else {
153  edm::LogPrint( "ParticleNotFound" )
154  << "BPHKinematicFit::kinParticles: "
155  << pname << " not found";
156  }
157  }
158  return plist;
159 }
160 
161 
163  if ( oldFit ) return kinematicTree( "", massConst, massSigma );
164  return kinTree;
165 }
166 
167 
169  const string& name,
170  double mass, double sigma ) const {
171  if ( sigma < 0 ) return kinematicTree( name, mass );
172  ParticleMass mc = mass;
173  MassKinematicConstraint kinConst( mc, sigma );
174  return kinematicTree( name, &kinConst );
175 }
176 
177 
179  const string& name,
180  double mass ) const {
181  if ( mass < 0 ) {
183  oldFit = false;
184  return kinTree;
185  }
186  int nn = daughFull().size();
187  ParticleMass mc = mass;
188  if ( nn == 2 ) {
189  TwoTrackMassKinematicConstraint kinConst( mc );
190  return kinematicTree( name, &kinConst );
191  }
192  else {
193  MultiTrackMassKinematicConstraint kinConst( mc, nn );
194  return kinematicTree( name, &kinConst );
195  }
196 }
197 
198 
200  const string& name,
201  KinematicConstraint* kc ) const {
203  oldFit = false;
204  kinParticles();
205  if ( allParticles.size() != daughFull().size() ) return kinTree;
206  vector<RefCountedKinematicParticle> kComp;
207  vector<RefCountedKinematicParticle> kTail;
208  if ( name != "" ) {
209  const BPHRecoCandidate* comp = getComp( name ).get();
210  if ( comp == 0 ) {
211  edm::LogPrint( "ParticleNotFound" )
212  << "BPHKinematicFit::kinematicTree: "
213  << name << " daughter not found";
214  return kinTree;
215  }
216  const vector<string>& names = comp->daugNames();
217  int ns;
218  int nn = ns = names.size();
219  vector<string> nfull( nn + 1 );
220  nfull[nn] = "*";
221  while ( nn-- ) nfull[nn] = name + "/" + names[nn];
222  vector<RefCountedKinematicParticle> kPart = kinParticles( nfull );
223  vector<RefCountedKinematicParticle>::const_iterator iter = kPart.begin();
224  vector<RefCountedKinematicParticle>::const_iterator imid = iter + ns;
225  vector<RefCountedKinematicParticle>::const_iterator iend = kPart.end();
226  kComp.insert( kComp.end(), iter, imid );
227  kTail.insert( kTail.end(), imid, iend );
228  }
229  else {
230  kComp = allParticles;
231  }
232  try {
234  RefCountedKinematicTree compTree = vtxFitter.fit( kComp );
235  if ( compTree->isEmpty() ) return kinTree;
236  KinematicParticleFitter kinFitter;
237  compTree = kinFitter.fit( kc, compTree );
238  if ( compTree->isEmpty() ) return kinTree;
239  compTree->movePointerToTheTop();
240  if ( kTail.size() ) {
241  RefCountedKinematicParticle compPart = compTree->currentParticle();
242  if ( !compPart->currentState().isValid() ) return kinTree;
243  kTail.push_back( compPart );
244  kinTree = vtxFitter.fit( kTail );
245  }
246  else {
247  kinTree = compTree;
248  }
249  }
250  catch ( std::exception e ) {
251  edm::LogPrint( "FitFailed" )
252  << "BPHKinematicFit::kinematicTree: "
253  << "kin fit reset";
255  }
256  return kinTree;
257 }
258 
259 
261  const string& name,
262  MultiTrackKinematicConstraint* kc ) const {
264  oldFit = false;
265  kinParticles();
266  if ( allParticles.size() != daughFull().size() ) return kinTree;
267  vector<string> nfull;
268  if ( name != "" ) {
269  const BPHRecoCandidate* comp = getComp( name ).get();
270  if ( comp == 0 ) {
271  edm::LogPrint( "ParticleNotFound" )
272  << "BPHKinematicFit::kinematicTree: "
273  << name << " daughter not found";
274  return kinTree;
275  }
276  const vector<string>& names = comp->daugNames();
277  int nn = names.size();
278  nfull.resize( nn + 1 );
279  nfull[nn] = "*";
280  while ( nn-- ) nfull[nn] = name + "/" + names[nn];
281  }
282  else {
283  nfull.push_back( "*" );
284  }
285  try {
287  kinTree = cvf.fit( kinParticles( nfull ), kc );
288  }
289  catch ( std::exception e ) {
290  edm::LogPrint( "FitFailed" )
291  << "BPHKinematicFit::kinematicTree: "
292  << "kin fit reset";
294  }
295  return kinTree;
296 }
297 
298 
300  oldKPs = oldFit = oldMom = true;
301  return;
302 }
303 
304 
306  kinematicTree();
307  if ( kinTree.get() == 0 ) return true;
308  return kinTree->isEmpty();
309 }
310 
311 
314  if ( kPart.get() == 0 ) return false;
315  return kPart->currentState().isValid();
316 }
317 
318 
320  if ( isEmpty() ) return RefCountedKinematicParticle( 0 );
321  return kinTree->currentParticle();
322 }
323 
324 
326  if ( isEmpty() ) return RefCountedKinematicVertex( 0 );
327  return kinTree->currentDecayVertex();
328 }
329 
330 
333  if ( kPart.get() == 0 ) return -1.0;
334  const KinematicState kStat = kPart->currentState();
335  if ( kStat.isValid() ) return kStat.mass();
336  return -1.0;
337 }
338 
339 
341  if ( oldMom ) fitMomentum();
342  return totalMomentum;
343 }
344 
345 
346 void BPHKinematicFit::addK( const string& name,
347  const reco::Candidate* daug,
348  double mass, double sigma ) {
349  addK( name, daug, "cfhpmig", mass, sigma );
350  return;
351 }
352 
353 
354 void BPHKinematicFit::addK( const string& name,
355  const reco::Candidate* daug,
356  const string& searchList,
357  double mass, double sigma ) {
358  addV( name, daug, searchList, mass );
359  dMSig[daughters().back()] = sigma;
360  return;
361 }
362 
363 
364 void BPHKinematicFit::addK( const string& name,
365  const BPHRecoConstCandPtr& comp ) {
366  addV( name, comp );
367  const map<const reco::Candidate*,double>& dMap = comp->dMSig;
368  dMSig.insert( dMap.begin(), dMap.end() );
369  return;
370 }
371 
372 
376  return;
377 }
378 
379 
381  kinMap.clear();
382  allParticles.clear();
383  const vector<const reco::Candidate*>& daug = daughFull();
385  int n = daug.size();
386  allParticles.reserve( n );
387  float chi = 0.0;
388  float ndf = 0.0;
389  while ( n-- ) {
390  const reco::Candidate* cand = daug[n];
391  ParticleMass mass = cand->mass();
392  float sigma = dMSig.find( cand )->second;
393  if ( sigma < 0 ) sigma = 1.0e-7;
395  if ( tt != 0 ) allParticles.push_back( kinMap[cand] =
396  pFactory.particle( *tt,
397  mass, chi, ndf, sigma ) );
398  }
399  oldKPs = false;
400  return;
401 }
402 
403 
405  if ( isValidFit() ) {
406  const KinematicState& ks = currentParticle()->currentState();
407  GlobalVector tm = ks.globalMomentum();
408  double x = tm.x();
409  double y = tm.y();
410  double z = tm.z();
411  double m = ks.mass();
412  double e = sqrt( ( x * x ) + ( y * y ) + ( z * z ) + ( m * m ) );
413  totalMomentum.SetPxPyPzE( x, y, z, e );
414  }
415  else {
416  edm::LogPrint( "FitNotFound" )
417  << "BPHKinematicFit::fitMomentum: "
418  << "simple momentum sum computed";
420  const vector<const reco::Candidate*>& daug = daughters();
421  int n = daug.size();
422  while ( n-- ) tm += daug[n]->p4();
423  const vector<BPHRecoConstCandPtr>& comp = daughComp();
424  int m = comp.size();
425  while ( m-- ) tm += comp[m]->p4();
426  totalMomentum = tm;
427  }
428  oldMom = false;
429  return;
430 }
431 
const reco::Candidate * cand
virtual const math::XYZTLorentzVector & p4() const
compute total momentum after the fit
RefCountedKinematicTree fit(const std::vector< RefCountedKinematicParticle > &part)
virtual ~BPHKinematicFit()
bool isValid() const
virtual void buildParticles() const
virtual void addK(const std::string &name, const reco::Candidate *daug, double mass=-1.0, double sigma=-1.0)
static const HistoName names[]
BPHGenericPtr< const BPHRecoCandidate >::type BPHRecoConstCandPtr
int kp
reco::TransientTrack * getTransientTrack(const reco::Candidate *cand) const
get TransientTrack for a daughter
double ParticleMass
Definition: ParticleMass.h:5
T y() const
Definition: PV3DBase.h:63
virtual bool isEmpty() const
GlobalVector globalMomentum() const
virtual ParticleMass mass() const
virtual const std::vector< const reco::Candidate * > & daughters() const
virtual const RefCountedKinematicParticle currentParticle() const
std::vector< RefCountedKinematicTree > fit(KinematicConstraint *cs, const std::vector< RefCountedKinematicTree > &trees) const
std::map< const reco::Candidate *, double > dMSig
ParticleMass mass() const
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
double constrSigma() const
std::map< const reco::Candidate *, RefCountedKinematicParticle > kinMap
double constrMass() const
retrieve the constraint
virtual void fitMomentum() const
virtual const RefCountedKinematicTree & kinematicTree() const
perform the kinematic fit and get the result
virtual void addV(const std::string &name, const reco::Candidate *daug, const std::string &searchList, double mass)
T sqrt(T t)
Definition: SSEVec.h:18
virtual const RefCountedKinematicVertex currentDecayVertex() const
virtual const std::vector< RefCountedKinematicParticle > & kinParticles() const
get kinematic particles
T z() const
Definition: PV3DBase.h:64
const std::vector< Component > & componentList() const
RefCountedKinematicTree fit(const std::vector< RefCountedKinematicParticle > &particles) const
virtual const reco::Candidate * originalReco(const reco::Candidate *daug) const
get the original particle from the clone
virtual void resetKinematicFit() const
reset the kinematic fit
ReferenceCountingPointer< KinematicVertex > RefCountedKinematicVertex
ReferenceCountingPointer< KinematicTree > RefCountedKinematicTree
std::vector< RefCountedKinematicParticle > allParticles
virtual bool isValidFit() const
RefCountedKinematicTree kinTree
virtual double mass() const =0
mass
std::map< std::string, const reco::Candidate * > dMap
ReferenceCountingPointer< KinematicParticle > RefCountedKinematicParticle
virtual BPHRecoConstCandPtr getComp(const std::string &name) const
virtual const std::vector< BPHRecoConstCandPtr > & daughComp() const
virtual void setNotUpdated() const
virtual const std::vector< const reco::Candidate * > & daughFull() const
RefCountedKinematicParticle particle(const reco::TransientTrack &initialTrack, const ParticleMass &massGuess, float chiSquared, float degreesOfFr, float &m_sigma) const
void setConstraint(double mass, double sigma)
apply a mass constraint
math::XYZTLorentzVector totalMomentum
virtual const std::vector< std::string > & daugNames() const
virtual const reco::Candidate * getDaug(const std::string &name) const
T x() const
Definition: PV3DBase.h:62
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger list("!*","!HLTx*"if it matches 2 triggers or more) will accept the event if all the matching triggers are FAIL.It will reject the event if any of the triggers are PASS or EXCEPTION(this matches the behavior of"!*"before the partial wildcard feature was incorporated).Triggers which are in the READY state are completely ignored.(READY should never be returned since the trigger paths have been run
virtual void setNotUpdated() const