CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
V0Fitter Class Reference

#include <RecoVertex/V0Producer/src/V0Fitter.cc>

Public Member Functions

const
reco::VertexCompositeCandidateCollection
getKshorts () const
 
const
reco::VertexCompositeCandidateCollection
getLambdas () const
 
 V0Fitter (const edm::ParameterSet &theParams, const edm::Event &iEvent, const edm::EventSetup &iSetup)
 
 ~V0Fitter ()
 

Private Member Functions

void cleanupFileOutput ()
 
double findV0MassError (const GlobalPoint &vtxPos, std::vector< reco::TransientTrack > dauTracks)
 
void fitAll (const edm::Event &iEvent, const edm::EventSetup &iSetup)
 
void initFileOutput ()
 

Private Attributes

double chi2Cut
 
double collinCut
 
bool doKshorts
 
bool doLambdas
 
double impactParameterSigCut
 
double innerHitPosCut
 
double kShortMassCut
 
double lambdaMassCut
 
const MagneticFieldmagField
 
double mPiPiCut
 
std::ofstream mPiPiMassOut
 
std::vector
< reco::TrackBase::TrackQuality
qualities
 
edm::InputTag recoAlg
 
double rVtxCut
 
bool storeRefTrax
 
reco::VertexCompositeCandidateCollection theKshorts
 
reco::VertexCompositeCandidateCollection theLambdas
 
double tkChi2Cut
 
double tkDCACut
 
int tkNhitsCut
 
const TrackerGeometrytrackerGeom
 
bool useRefTrax
 
edm::InputTag vtxFitter
 
double vtxSigCut
 

Detailed Description

Description: <one line="" class="" summary>="">

Implementation: <Notes on="" implementation>="">

Definition at line 57 of file V0Fitter.h.

Constructor & Destructor Documentation

V0Fitter::V0Fitter ( const edm::ParameterSet theParams,
const edm::Event iEvent,
const edm::EventSetup iSetup 
)

Definition at line 48 of file V0Fitter.cc.

References chi2Cut, collinCut, doKshorts, doLambdas, fitAll(), edm::ParameterSet::getParameter(), impactParameterSigCut, innerHitPosCut, kShortMassCut, lambdaMassCut, mPiPiCut, qualities, reco::TrackBase::qualityByName(), recoAlg, rVtxCut, tkChi2Cut, tkDCACut, tkNhitsCut, useRefTrax, vtxFitter, and vtxSigCut.

49  {
50  using std::string;
51 
52  // Get the track reco algorithm from the ParameterSet
53  recoAlg = theParameters.getParameter<edm::InputTag>("trackRecoAlgorithm");
54 
55  // ------> Initialize parameters from PSet. ALL TRACKED, so no defaults.
56  // First set bits to do various things:
57  // -decide whether to use the KVF track smoother, and whether to store those
58  // tracks in the reco::Vertex
59  useRefTrax = theParameters.getParameter<bool>(string("useSmoothing"));
60 
61  // -whether to reconstruct K0s
62  doKshorts = theParameters.getParameter<bool>(string("selectKshorts"));
63  // -whether to reconstruct Lambdas
64  doLambdas = theParameters.getParameter<bool>(string("selectLambdas"));
65 
66  // Second, initialize post-fit cuts
67  chi2Cut = theParameters.getParameter<double>(string("vtxChi2Cut"));
68  tkChi2Cut = theParameters.getParameter<double>(string("tkChi2Cut"));
69  tkNhitsCut = theParameters.getParameter<int>(string("tkNhitsCut"));
70  rVtxCut = theParameters.getParameter<double>(string("rVtxCut"));
71  vtxSigCut = theParameters.getParameter<double>(string("vtxSignificance2DCut"));
72  collinCut = theParameters.getParameter<double>(string("collinearityCut"));
73  kShortMassCut = theParameters.getParameter<double>(string("kShortMassCut"));
74  lambdaMassCut = theParameters.getParameter<double>(string("lambdaMassCut"));
75  impactParameterSigCut = theParameters.getParameter<double>(string("impactParameterSigCut"));
76  mPiPiCut = theParameters.getParameter<double>(string("mPiPiCut"));
77  tkDCACut = theParameters.getParameter<double>(string("tkDCACut"));
78  vtxFitter = theParameters.getParameter<edm::InputTag>("vertexFitter");
79  innerHitPosCut = theParameters.getParameter<double>(string("innerHitPosCut"));
80  std::vector<std::string> qual = theParameters.getParameter<std::vector<std::string> >("trackQualities");
81  for (unsigned int ndx = 0; ndx < qual.size(); ndx++) {
82  qualities.push_back(reco::TrackBase::qualityByName(qual[ndx]));
83  }
84 
85  //edm::LogInfo("V0Producer") << "Using " << vtxFitter << " to fit V0 vertices.\n";
86  //std::cout << "Using " << vtxFitter << " to fit V0 vertices." << std::endl;
87  // FOR DEBUG:
88  //initFileOutput();
89  //--------------------
90 
91  //std::cout << "Entering V0Producer" << std::endl;
92 
93  fitAll(iEvent, iSetup);
94 
95  // FOR DEBUG:
96  //cleanupFileOutput();
97  //--------------------
98 
99 }
bool useRefTrax
Definition: V0Fitter.h:78
int tkNhitsCut
Definition: V0Fitter.h:89
void fitAll(const edm::Event &iEvent, const edm::EventSetup &iSetup)
Definition: V0Fitter.cc:105
double chi2Cut
Definition: V0Fitter.h:87
edm::InputTag vtxFitter
Definition: V0Fitter.h:102
double rVtxCut
Definition: V0Fitter.h:90
double collinCut
Definition: V0Fitter.h:92
bool doKshorts
Definition: V0Fitter.h:80
double tkChi2Cut
Definition: V0Fitter.h:88
double impactParameterSigCut
Definition: V0Fitter.h:95
bool doLambdas
Definition: V0Fitter.h:81
std::vector< reco::TrackBase::TrackQuality > qualities
Definition: V0Fitter.h:100
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:46
double tkDCACut
Definition: V0Fitter.h:97
double kShortMassCut
Definition: V0Fitter.h:93
double mPiPiCut
Definition: V0Fitter.h:96
edm::InputTag recoAlg
Definition: V0Fitter.h:77
double vtxSigCut
Definition: V0Fitter.h:91
double innerHitPosCut
Definition: V0Fitter.h:98
double lambdaMassCut
Definition: V0Fitter.h:94
V0Fitter::~V0Fitter ( )

Definition at line 101 of file V0Fitter.cc.

101  {
102 }

Member Function Documentation

void V0Fitter::cleanupFileOutput ( )
inlineprivate

Definition at line 117 of file V0Fitter.h.

References mPiPiMassOut.

117  {
118  mPiPiMassOut.close();
119  }
std::ofstream mPiPiMassOut
Definition: V0Fitter.h:112
double V0Fitter::findV0MassError ( const GlobalPoint vtxPos,
std::vector< reco::TransientTrack dauTracks 
)
private

Definition at line 513 of file V0Fitter.cc.

513  {
514  return -1.;
515 }
void V0Fitter::fitAll ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
private

Definition at line 105 of file V0Fitter.cc.

References abs, ClosestApproachInRPhi::calculate(), DeDxDiscriminatorTools::charge(), reco::Vertex::chi2(), chi2Cut, gather_cfg::cout, reco::Vertex::covariance(), ClosestApproachInRPhi::crossingPoint(), ClosestApproachInRPhi::distance(), doKshorts, doLambdas, edm::EventSetup::get(), edm::Event::getByLabel(), TransientVertex::hasRefittedTracks(), impactParameterSigCut, reco::TransientTrack::impactPointTSCP(), TrajectoryStateTransform::initialFreeState(), innerHitPosCut, TrajectoryStateClosestToBeamLine::isValid(), TransientVertex::isValid(), TrajectoryStateClosestToPoint::isValid(), kShortMass, kShortMassCut, lambdaMass, lambdaMassCut, PV3DBase< T, PVType, FrameType >::mag2(), mag2(), magField, TrajectoryStateClosestToPoint::momentum(), mPiPiCut, reco::Vertex::ndof(), reco::Vertex::normalizedChi2(), piMassSquared, protonMassSquared, qualities, dt_dqm_sourceclient_common_cff::reco, recoAlg, TransientVertex::refittedTracks(), rVtxCut, AddFourMomenta::set(), reco::RecoChargedCandidate::setTrack(), Measurement1D::significance(), mathSSE::sqrt(), ClosestApproachInRPhi::status(), theKshorts, theLambdas, TrajectoryStateClosestToPoint::theState(), tkChi2Cut, tkDCACut, tkNhitsCut, TransientVertex::totalChiSquared(), trackerGeom, reco::TransientTrack::trajectoryStateClosestToPoint(), TrajectoryStateClosestToBeamLine::transverseImpactParameter(), useRefTrax, KalmanVertexFitter::vertex(), AdaptiveVertexFitter::vertex(), vtxFitter, vtxSigCut, reco::Vertex::x(), reco::Vertex::y(), and reco::Vertex::z().

Referenced by V0Fitter().

105  {
106 
107  using std::vector;
108  using std::cout;
109  using std::endl;
110  using namespace reco;
111  using namespace edm;
112 
113  // Create std::vectors for Tracks and TrackRefs (required for
114  // passing to the KalmanVertexFitter)
115  std::vector<TrackRef> theTrackRefs;
116  std::vector<TransientTrack> theTransTracks;
117 
118  // Handles for tracks, B-field, and tracker geometry
119  Handle<reco::TrackCollection> theTrackHandle;
120  Handle<reco::BeamSpot> theBeamSpotHandle;
121  ESHandle<MagneticField> bFieldHandle;
122  ESHandle<TrackerGeometry> trackerGeomHandle;
123  ESHandle<GlobalTrackingGeometry> globTkGeomHandle;
124  //cout << "Check 0" << endl;
125 
126  // Get the tracks from the event, and get the B-field record
127  // from the EventSetup
128  iEvent.getByLabel(recoAlg, theTrackHandle);
129  iEvent.getByLabel(std::string("offlineBeamSpot"), theBeamSpotHandle);
130  if( !theTrackHandle->size() ) return;
131  iSetup.get<IdealMagneticFieldRecord>().get(bFieldHandle);
132  iSetup.get<TrackerDigiGeometryRecord>().get(trackerGeomHandle);
133  iSetup.get<GlobalTrackingGeometryRecord>().get(globTkGeomHandle);
134 
135  trackerGeom = trackerGeomHandle.product();
136  magField = bFieldHandle.product();
137 
138  // Fill vectors of TransientTracks and TrackRefs after applying preselection cuts.
139  for(unsigned int indx = 0; indx < theTrackHandle->size(); indx++) {
140  TrackRef tmpRef( theTrackHandle, indx );
141  bool quality_ok = true;
142  if (qualities.size()!=0) {
143  quality_ok = false;
144  for (unsigned int ndx_ = 0; ndx_ < qualities.size(); ndx_++) {
145  if (tmpRef->quality(qualities[ndx_])){
146  quality_ok = true;
147  break;
148  }
149  }
150  }
151  if( !quality_ok ) continue;
152 
153 
154  if( tmpRef->normalizedChi2() < tkChi2Cut &&
155  tmpRef->numberOfValidHits() >= tkNhitsCut ) {
156  TransientTrack tmpTk( *tmpRef, &(*bFieldHandle), globTkGeomHandle );
157  TrajectoryStateTransform theTransform;
158  FreeTrajectoryState initialFTS = theTransform.initialFreeState(*tmpRef, magField);
159  TSCBLBuilderNoMaterial blsBuilder;
160  TrajectoryStateClosestToBeamLine tscb( blsBuilder(initialFTS, *theBeamSpotHandle) );
161 
162  if( tscb.isValid() ) {
163  if( tscb.transverseImpactParameter().significance() > impactParameterSigCut ) {
164  theTrackRefs.push_back( tmpRef );
165  theTransTracks.push_back( tmpTk );
166  }
167  }
168  }
169  }
170 
171  // Good tracks have now been selected for vertexing. Move on to vertex fitting.
172 
173 
174  // Loop over tracks and vertex good charged track pairs
175  for(unsigned int trdx1 = 0; trdx1 < theTrackRefs.size(); trdx1++) {
176 
177  for(unsigned int trdx2 = trdx1 + 1; trdx2 < theTrackRefs.size(); trdx2++) {
178 
179  //This vector holds the pair of oppositely-charged tracks to be vertexed
180  std::vector<TransientTrack> transTracks;
181 
182  TrackRef positiveTrackRef;
183  TrackRef negativeTrackRef;
184  TransientTrack* posTransTkPtr = 0;
185  TransientTrack* negTransTkPtr = 0;
186 
187  // Look at the two tracks we're looping over. If they're oppositely
188  // charged, load them into the hypothesized positive and negative tracks
189  // and references to be sent to the KalmanVertexFitter
190  if(theTrackRefs[trdx1]->charge() < 0. &&
191  theTrackRefs[trdx2]->charge() > 0.) {
192  negativeTrackRef = theTrackRefs[trdx1];
193  positiveTrackRef = theTrackRefs[trdx2];
194  negTransTkPtr = &theTransTracks[trdx1];
195  posTransTkPtr = &theTransTracks[trdx2];
196  }
197  else if(theTrackRefs[trdx1]->charge() > 0. &&
198  theTrackRefs[trdx2]->charge() < 0.) {
199  negativeTrackRef = theTrackRefs[trdx2];
200  positiveTrackRef = theTrackRefs[trdx1];
201  negTransTkPtr = &theTransTracks[trdx2];
202  posTransTkPtr = &theTransTracks[trdx1];
203  }
204  // If they're not 2 oppositely charged tracks, loop back to the
205  // beginning and try the next pair.
206  else continue;
207 
208  // Fill the vector of TransientTracks to send to KVF
209  transTracks.push_back(*posTransTkPtr);
210  transTracks.push_back(*negTransTkPtr);
211 
212  // Trajectory states to calculate DCA for the 2 tracks
213  FreeTrajectoryState posState = posTransTkPtr->impactPointTSCP().theState();
214  FreeTrajectoryState negState = negTransTkPtr->impactPointTSCP().theState();
215 
216  if( !posTransTkPtr->impactPointTSCP().isValid() || !negTransTkPtr->impactPointTSCP().isValid() ) continue;
217 
218  // Measure distance between tracks at their closest approach
220  cApp.calculate(posState, negState);
221  if( !cApp.status() ) continue;
222  float dca = fabs( cApp.distance() );
223  GlobalPoint cxPt = cApp.crossingPoint();
224 
225  if (dca < 0. || dca > tkDCACut) continue;
226  if (sqrt( cxPt.x()*cxPt.x() + cxPt.y()*cxPt.y() ) > 120.
227  || std::abs(cxPt.z()) > 300.) continue;
228 
229  // Get trajectory states for the tracks at POCA for later cuts
231  posTransTkPtr->trajectoryStateClosestToPoint( cxPt );
233  negTransTkPtr->trajectoryStateClosestToPoint( cxPt );
234 
235  if( !posTSCP.isValid() || !negTSCP.isValid() ) continue;
236 
237 
238  /*double posESq = posTSCP.momentum().mag2() + piMassSquared;
239  double negESq = negTSCP.momentum().mag2() + piMassSquared;
240  double posE = sqrt(posESq);
241  double negE = sqrt(negESq);
242  double totalE = posE + negE;*/
243  double totalE = sqrt( posTSCP.momentum().mag2() + piMassSquared ) +
244  sqrt( negTSCP.momentum().mag2() + piMassSquared );
245  double totalESq = totalE*totalE;
246  double totalPSq =
247  ( posTSCP.momentum() + negTSCP.momentum() ).mag2();
248  double mass = sqrt( totalESq - totalPSq);
249 
250  //mPiPiMassOut << mass << std::endl;
251 
252  if( mass > mPiPiCut ) continue;
253 
254  // Create the vertex fitter object and vertex the tracks
255  TransientVertex theRecoVertex;
256  if(vtxFitter == std::string("KalmanVertexFitter")) {
257  KalmanVertexFitter theKalmanFitter(useRefTrax == 0 ? false : true);
258  theRecoVertex = theKalmanFitter.vertex(transTracks);
259  }
260  else if (vtxFitter == std::string("AdaptiveVertexFitter")) {
261  useRefTrax = false;
262  AdaptiveVertexFitter theAdaptiveFitter;
263  theRecoVertex = theAdaptiveFitter.vertex(transTracks);
264  }
265 
266  // If the vertex is valid, make a VertexCompositeCandidate with it
267 
268  if( !theRecoVertex.isValid() || theRecoVertex.totalChiSquared() < 0. ) {
269  continue;
270  }
271 
272  // Create reco::Vertex object for use in creating the Candidate
273  reco::Vertex theVtx = theRecoVertex;
274  // Create and fill vector of refitted TransientTracks
275  // (iff they've been created by the KVF)
276  std::vector<TransientTrack> refittedTrax;
277  if( theRecoVertex.hasRefittedTracks() ) {
278  refittedTrax = theRecoVertex.refittedTracks();
279  }
280 
281  // Do post-fit cuts if specified in config file.
282 
283  // Find the vertex d0 and its error
284 
285  typedef ROOT::Math::SMatrix<double, 3, 3, ROOT::Math::MatRepSym<double, 3> > SMatrixSym3D;
286  typedef ROOT::Math::SVector<double, 3> SVector3;
287 
288  GlobalPoint vtxPos(theVtx.x(), theVtx.y(), theVtx.z());
289 
290  GlobalPoint beamSpotPos(theBeamSpotHandle->position().x(),
291  theBeamSpotHandle->position().y(),
292  theBeamSpotHandle->position().z());
293 
294  SMatrixSym3D totalCov = theBeamSpotHandle->covariance3D() + theVtx.covariance();
295  SVector3 distanceVector(vtxPos.x() - beamSpotPos.x(),
296  vtxPos.y() - beamSpotPos.y(),
297  0.);//so that we get radial values only,
298  //since z beamSpot uncertainty is huge
299 
300  double rVtxMag = ROOT::Math::Mag(distanceVector);
301  double sigmaRvtxMag = sqrt(ROOT::Math::Similarity(totalCov, distanceVector)) / rVtxMag;
302 
303  // The methods innerOk() and innerPosition() require TrackExtra, which
304  // is only available in the RECO data tier, not AOD. Setting innerHitPosCut
305  // to -1 avoids this problem and allows to run on AOD.
306  if( innerHitPosCut > 0. && positiveTrackRef->innerOk() ) {
307  reco::Vertex::Point posTkHitPos = positiveTrackRef->innerPosition();
308  double posTkHitPosD2 =
309  (posTkHitPos.x()-beamSpotPos.x())*(posTkHitPos.x()-beamSpotPos.x()) +
310  (posTkHitPos.y()-beamSpotPos.y())*(posTkHitPos.y()-beamSpotPos.y());
311  if( sqrt( posTkHitPosD2 ) < ( rVtxMag - sigmaRvtxMag*innerHitPosCut )
312  ) {
313  continue;
314  }
315  }
316  if( innerHitPosCut > 0. && negativeTrackRef->innerOk() ) {
317  reco::Vertex::Point negTkHitPos = negativeTrackRef->innerPosition();
318  double negTkHitPosD2 =
319  (negTkHitPos.x()-beamSpotPos.x())*(negTkHitPos.x()-beamSpotPos.x()) +
320  (negTkHitPos.y()-beamSpotPos.y())*(negTkHitPos.y()-beamSpotPos.y());
321  if( sqrt( negTkHitPosD2 ) < ( rVtxMag - sigmaRvtxMag*innerHitPosCut )
322  ) {
323  continue;
324  }
325  }
326 
327  if( theVtx.normalizedChi2() > chi2Cut ||
328  rVtxMag < rVtxCut ||
329  rVtxMag / sigmaRvtxMag < vtxSigCut ) {
330  continue;
331  }
332 
333  // Cuts finished, now we create the candidates and push them back into the collections.
334 
337 
338  if( useRefTrax && refittedTrax.size() > 1 ) {
339  // Need an iterator over the refitted tracks for below
340  std::vector<TransientTrack>::iterator traxIter = refittedTrax.begin(),
341  traxEnd = refittedTrax.end();
342 
343  // TransientTrack objects to hold the positive and negative
344  // refitted tracks
345  TransientTrack* thePositiveRefTrack = 0;
346  TransientTrack* theNegativeRefTrack = 0;
347 
348  for( ; traxIter != traxEnd; ++traxIter) {
349  if( traxIter->track().charge() > 0. ) {
350  thePositiveRefTrack = new TransientTrack(*traxIter);
351  }
352  else if (traxIter->track().charge() < 0.) {
353  theNegativeRefTrack = new TransientTrack(*traxIter);
354  }
355  }
356  trajPlus = new TrajectoryStateClosestToPoint(
357  thePositiveRefTrack->trajectoryStateClosestToPoint(vtxPos));
358  trajMins = new TrajectoryStateClosestToPoint(
359  theNegativeRefTrack->trajectoryStateClosestToPoint(vtxPos));
360  delete thePositiveRefTrack;
361  delete theNegativeRefTrack;
362  }
363  else {
364  trajPlus = new TrajectoryStateClosestToPoint(
365  posTransTkPtr->trajectoryStateClosestToPoint(vtxPos));
366  trajMins = new TrajectoryStateClosestToPoint(
367  negTransTkPtr->trajectoryStateClosestToPoint(vtxPos));
368 
369  }
370 
371  if( !trajPlus->isValid() || !trajMins->isValid() ) continue;
372 
373  posTransTkPtr = negTransTkPtr = 0;
374 
375  GlobalVector positiveP(trajPlus->momentum());
376  GlobalVector negativeP(trajMins->momentum());
377  GlobalVector totalP(positiveP + negativeP);
378 
379  //cleanup stuff we don't need anymore
380  delete trajPlus;
381  delete trajMins;
382  trajPlus = trajMins = 0;
383 
384  // calculate total energy of V0 3 ways:
385  // Assume it's a kShort, a Lambda, or a LambdaBar.
386  double piPlusE = sqrt( positiveP.mag2() + piMassSquared );
387  double piMinusE = sqrt( negativeP.mag2() + piMassSquared );
388  double protonE = sqrt( positiveP.mag2() + protonMassSquared );
389  double antiProtonE = sqrt( negativeP.mag2() + protonMassSquared );
390  double kShortETot = piPlusE + piMinusE;
391  double lambdaEtot = protonE + piMinusE;
392  double lambdaBarEtot = antiProtonE + piPlusE;
393 
394  using namespace reco;
395 
396  // Create momentum 4-vectors for the 3 candidate types
397  const Particle::LorentzVector kShortP4(totalP.x(),
398  totalP.y(), totalP.z(),
399  kShortETot);
400  const Particle::LorentzVector lambdaP4(totalP.x(),
401  totalP.y(), totalP.z(),
402  lambdaEtot);
403  const Particle::LorentzVector lambdaBarP4(totalP.x(),
404  totalP.y(), totalP.z(),
405  lambdaBarEtot);
406 
407  Particle::Point vtx(theVtx.x(), theVtx.y(), theVtx.z());
408  const Vertex::CovarianceMatrix vtxCov(theVtx.covariance());
409  double vtxChi2(theVtx.chi2());
410  double vtxNdof(theVtx.ndof());
411 
412  // Create the VertexCompositeCandidate object that will be stored in the Event
413  VertexCompositeCandidate* theKshort = 0;
414  VertexCompositeCandidate* theLambda = 0;
415  VertexCompositeCandidate* theLambdaBar = 0;
416 
417  if( doKshorts ) {
418  theKshort = new VertexCompositeCandidate(0, kShortP4, vtx, vtxCov, vtxChi2, vtxNdof);
419  }
420  if( doLambdas ) {
421  if( positiveP.mag() > negativeP.mag() ) {
422  theLambda =
423  new VertexCompositeCandidate(0, lambdaP4, vtx, vtxCov, vtxChi2, vtxNdof);
424  }
425  else {
426  theLambdaBar =
427  new VertexCompositeCandidate(0, lambdaBarP4, vtx, vtxCov, vtxChi2, vtxNdof);
428  }
429  }
430 
431  // Create daughter candidates for the VertexCompositeCandidates
433  thePiPlusCand(1, Particle::LorentzVector(positiveP.x(),
434  positiveP.y(), positiveP.z(),
435  piPlusE), vtx);
436  thePiPlusCand.setTrack(positiveTrackRef);
437 
439  thePiMinusCand(-1, Particle::LorentzVector(negativeP.x(),
440  negativeP.y(), negativeP.z(),
441  piMinusE), vtx);
442  thePiMinusCand.setTrack(negativeTrackRef);
443 
444 
446  theProtonCand(1, Particle::LorentzVector(positiveP.x(),
447  positiveP.y(), positiveP.z(),
448  protonE), vtx);
449  theProtonCand.setTrack(positiveTrackRef);
450 
452  theAntiProtonCand(-1, Particle::LorentzVector(negativeP.x(),
453  negativeP.y(), negativeP.z(),
454  antiProtonE), vtx);
455  theAntiProtonCand.setTrack(negativeTrackRef);
456 
457 
458  AddFourMomenta addp4;
459  // Store the daughter Candidates in the VertexCompositeCandidates
460  // if they pass mass cuts
461  if( doKshorts ) {
462  theKshort->addDaughter(thePiPlusCand);
463  theKshort->addDaughter(thePiMinusCand);
464  theKshort->setPdgId(310);
465  addp4.set( *theKshort );
466  if( theKshort->mass() < kShortMass + kShortMassCut &&
467  theKshort->mass() > kShortMass - kShortMassCut ) {
468  theKshorts.push_back( *theKshort );
469  }
470  }
471 
472  if( doLambdas && theLambda ) {
473  theLambda->addDaughter(theProtonCand);
474  theLambda->addDaughter(thePiMinusCand);
475  theLambda->setPdgId(3122);
476  addp4.set( *theLambda );
477  if( theLambda->mass() < lambdaMass + lambdaMassCut &&
478  theLambda->mass() > lambdaMass - lambdaMassCut ) {
479  theLambdas.push_back( *theLambda );
480  }
481  }
482  else if ( doLambdas && theLambdaBar ) {
483  theLambdaBar->addDaughter(theAntiProtonCand);
484  theLambdaBar->addDaughter(thePiPlusCand);
485  theLambdaBar->setPdgId(-3122);
486  addp4.set( *theLambdaBar );
487  if( theLambdaBar->mass() < lambdaMass + lambdaMassCut &&
488  theLambdaBar->mass() > lambdaMass - lambdaMassCut ) {
489  theLambdas.push_back( *theLambdaBar );
490  }
491  }
492 
493  if(theKshort) delete theKshort;
494  if(theLambda) delete theLambda;
495  if(theLambdaBar) delete theLambdaBar;
496  theKshort = theLambda = theLambdaBar = 0;
497 
498  }
499  }
500 }
bool useRefTrax
Definition: V0Fitter.h:78
T mag2() const
Definition: PV3DBase.h:60
reco::VertexCompositeCandidateCollection theLambdas
Definition: V0Fitter.h:70
const double lambdaMass
Definition: V0Fitter.cc:45
int tkNhitsCut
Definition: V0Fitter.h:89
TrajectoryStateClosestToPoint impactPointTSCP() const
const FreeTrajectoryState & theState() const
double y() const
y coordinate
Definition: Vertex.h:97
float totalChiSquared() const
const MagneticField * magField
Definition: V0Fitter.h:75
double chi2Cut
Definition: V0Fitter.h:87
#define abs(x)
Definition: mlp_lapack.h:159
const double piMassSquared
Definition: V0Fitter.cc:41
const double kShortMass
Definition: V0Fitter.cc:44
double covariance(int i, int j) const
(i, j)-th element of error matrix, i, j = 0, ... 2
Definition: Vertex.h:110
bool hasRefittedTracks() const
edm::InputTag vtxFitter
Definition: V0Fitter.h:102
double rVtxCut
Definition: V0Fitter.h:90
double charge(const std::vector< uint8_t > &Ampls)
math::Error< dimension >::type CovarianceMatrix
covariance error matrix (3x3)
Definition: Vertex.h:46
virtual GlobalPoint crossingPoint() const
bool doKshorts
Definition: V0Fitter.h:80
double tkChi2Cut
Definition: V0Fitter.h:88
double impactParameterSigCut
Definition: V0Fitter.h:95
T sqrt(T t)
Definition: SSEVec.h:28
bool doLambdas
Definition: V0Fitter.h:81
const double protonMassSquared
Definition: V0Fitter.cc:43
T mag2() const
The vector magnitude squared. Equivalent to vec.dot(vec)
const TrackerGeometry * trackerGeom
Definition: V0Fitter.h:73
virtual CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &) const
math::XYZPoint Point
point in the space
Definition: Particle.h:30
double chi2() const
chi-squares
Definition: Vertex.h:82
double z() const
y coordinate
Definition: Vertex.h:99
math::XYZPoint Point
point in the space
Definition: Vertex.h:40
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:355
std::vector< reco::TransientTrack > refittedTracks() const
std::vector< reco::TrackBase::TrackQuality > qualities
Definition: V0Fitter.h:100
double ndof() const
Definition: Vertex.h:89
double x() const
x coordinate
Definition: Vertex.h:95
virtual bool calculate(const TrajectoryStateOnSurface &sta, const TrajectoryStateOnSurface &stb)
double tkDCACut
Definition: V0Fitter.h:97
TrajectoryStateClosestToPoint trajectoryStateClosestToPoint(const GlobalPoint &point) const
const T & get() const
Definition: EventSetup.h:55
reco::VertexCompositeCandidateCollection theKshorts
Definition: V0Fitter.h:69
double kShortMassCut
Definition: V0Fitter.h:93
double mPiPiCut
Definition: V0Fitter.h:96
edm::InputTag recoAlg
Definition: V0Fitter.h:77
tuple cout
Definition: gather_cfg.py:41
void set(reco::Candidate &c) const
set up a candidate
double normalizedChi2() const
chi-squared divided by n.d.o.f.
Definition: Vertex.h:91
double vtxSigCut
Definition: V0Fitter.h:91
virtual float distance() const
double innerHitPosCut
Definition: V0Fitter.h:98
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:26
bool isValid() const
virtual bool status() const
double lambdaMassCut
Definition: V0Fitter.h:94
FreeTrajectoryState initialFreeState(const reco::Track &tk, const MagneticField *field) const
const reco::VertexCompositeCandidateCollection & V0Fitter::getKshorts ( ) const

Definition at line 503 of file V0Fitter.cc.

References theKshorts.

Referenced by V0Producer::produce().

503  {
504  return theKshorts;
505 }
reco::VertexCompositeCandidateCollection theKshorts
Definition: V0Fitter.h:69
const reco::VertexCompositeCandidateCollection & V0Fitter::getLambdas ( ) const

Definition at line 507 of file V0Fitter.cc.

References theLambdas.

Referenced by V0Producer::produce().

507  {
508  return theLambdas;
509 }
reco::VertexCompositeCandidateCollection theLambdas
Definition: V0Fitter.h:70
void V0Fitter::initFileOutput ( )
inlineprivate

Definition at line 114 of file V0Fitter.h.

References mPiPiMassOut.

114  {
115  mPiPiMassOut.open("mPiPi.txt", std::ios::app);
116  }
std::ofstream mPiPiMassOut
Definition: V0Fitter.h:112

Member Data Documentation

double V0Fitter::chi2Cut
private

Definition at line 87 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::collinCut
private

Definition at line 92 of file V0Fitter.h.

Referenced by V0Fitter().

bool V0Fitter::doKshorts
private

Definition at line 80 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

bool V0Fitter::doLambdas
private

Definition at line 81 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::impactParameterSigCut
private

Definition at line 95 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::innerHitPosCut
private

Definition at line 98 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::kShortMassCut
private

Definition at line 93 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::lambdaMassCut
private

Definition at line 94 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

const MagneticField* V0Fitter::magField
private

Definition at line 75 of file V0Fitter.h.

Referenced by fitAll().

double V0Fitter::mPiPiCut
private

Definition at line 96 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

std::ofstream V0Fitter::mPiPiMassOut
private

Definition at line 112 of file V0Fitter.h.

Referenced by cleanupFileOutput(), and initFileOutput().

std::vector<reco::TrackBase::TrackQuality> V0Fitter::qualities
private

Definition at line 100 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

edm::InputTag V0Fitter::recoAlg
private

Definition at line 77 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::rVtxCut
private

Definition at line 90 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

bool V0Fitter::storeRefTrax
private

Definition at line 79 of file V0Fitter.h.

reco::VertexCompositeCandidateCollection V0Fitter::theKshorts
private

Definition at line 69 of file V0Fitter.h.

Referenced by fitAll(), and getKshorts().

reco::VertexCompositeCandidateCollection V0Fitter::theLambdas
private

Definition at line 70 of file V0Fitter.h.

Referenced by fitAll(), and getLambdas().

double V0Fitter::tkChi2Cut
private

Definition at line 88 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::tkDCACut
private

Definition at line 97 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

int V0Fitter::tkNhitsCut
private

Definition at line 89 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

const TrackerGeometry* V0Fitter::trackerGeom
private

Definition at line 73 of file V0Fitter.h.

Referenced by fitAll().

bool V0Fitter::useRefTrax
private

Definition at line 78 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

edm::InputTag V0Fitter::vtxFitter
private

Definition at line 102 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::vtxSigCut
private

Definition at line 91 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().