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

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

Private Member Functions

void cleanupFileOutput ()
 
double findV0MassError (const GlobalPoint &vtxPos, const std::vector< reco::TransientTrack > &dauTracks)
 
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
 
double rVtxCut
 
reco::VertexCompositeCandidateCollection theKshorts
 
reco::VertexCompositeCandidateCollection theLambdas
 
double tkChi2Cut
 
double tkDCACut
 
int tkNhitsCut
 
edm::EDGetTokenT< reco::BeamSpottoken_beamSpot
 
edm::EDGetTokenT
< reco::TrackCollection
token_tracks
 
const TrackerGeometrytrackerGeom
 
bool useRefTrax
 
edm::InputTag vtxFitter
 
double vtxSigCut
 

Detailed Description

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

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

Definition at line 58 of file V0Fitter.h.

Constructor & Destructor Documentation

V0Fitter::V0Fitter ( const edm::ParameterSet theParams,
edm::ConsumesCollector &&  iC 
)

Definition at line 46 of file V0Fitter.cc.

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

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

Definition at line 99 of file V0Fitter.cc.

99  {
100 }

Member Function Documentation

void V0Fitter::cleanupFileOutput ( )
inlineprivate

Definition at line 120 of file V0Fitter.h.

References mPiPiMassOut.

120  {
121  mPiPiMassOut.close();
122  }
std::ofstream mPiPiMassOut
Definition: V0Fitter.h:115
double V0Fitter::findV0MassError ( const GlobalPoint vtxPos,
const std::vector< reco::TransientTrack > &  dauTracks 
)
private

Definition at line 508 of file V0Fitter.cc.

508  {
509  return -1.;
510 }
void V0Fitter::fitAll ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)

Definition at line 103 of file V0Fitter.cc.

References funct::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::getByToken(), TransientVertex::hasRefittedTracks(), impactParameterSigCut, reco::TransientTrack::impactPointTSCP(), trajectoryStateTransform::initialFreeState(), innerHitPosCut, TrajectoryStateClosestToBeamLine::isValid(), TrajectoryStateClosestToPoint::isValid(), TransientVertex::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, TransientVertex::refittedTracks(), rVtxCut, AddFourMomenta::set(), reco::RecoChargedCandidate::setTrack(), Measurement1D::significance(), mathSSE::sqrt(), ClosestApproachInRPhi::status(), AlCaHLTBitMon_QueryRunRegistry::string, theKshorts, theLambdas, TrajectoryStateClosestToPoint::theState(), tkChi2Cut, tkDCACut, tkNhitsCut, token_beamSpot, token_tracks, 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 V0Producer::produce().

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

Definition at line 498 of file V0Fitter.cc.

References theKshorts.

Referenced by V0Producer::produce().

498  {
499  return theKshorts;
500 }
reco::VertexCompositeCandidateCollection theKshorts
Definition: V0Fitter.h:71
const reco::VertexCompositeCandidateCollection & V0Fitter::getLambdas ( ) const

Definition at line 502 of file V0Fitter.cc.

References theLambdas.

Referenced by V0Producer::produce().

502  {
503  return theLambdas;
504 }
reco::VertexCompositeCandidateCollection theLambdas
Definition: V0Fitter.h:72
void V0Fitter::initFileOutput ( )
inlineprivate

Definition at line 117 of file V0Fitter.h.

References mPiPiMassOut.

117  {
118  mPiPiMassOut.open("mPiPi.txt", std::ios::app);
119  }
std::ofstream mPiPiMassOut
Definition: V0Fitter.h:115

Member Data Documentation

double V0Fitter::chi2Cut
private

Definition at line 88 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::collinCut
private

Definition at line 94 of file V0Fitter.h.

Referenced by V0Fitter().

bool V0Fitter::doKshorts
private

Definition at line 81 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

bool V0Fitter::doLambdas
private

Definition at line 82 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::impactParameterSigCut
private

Definition at line 97 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::innerHitPosCut
private

Definition at line 100 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::kShortMassCut
private

Definition at line 95 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::lambdaMassCut
private

Definition at line 96 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

const MagneticField* V0Fitter::magField
private

Definition at line 77 of file V0Fitter.h.

Referenced by fitAll().

double V0Fitter::mPiPiCut
private

Definition at line 98 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

std::ofstream V0Fitter::mPiPiMassOut
private

Definition at line 115 of file V0Fitter.h.

Referenced by cleanupFileOutput(), and initFileOutput().

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

Definition at line 102 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::rVtxCut
private

Definition at line 91 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

reco::VertexCompositeCandidateCollection V0Fitter::theKshorts
private

Definition at line 71 of file V0Fitter.h.

Referenced by fitAll(), and getKshorts().

reco::VertexCompositeCandidateCollection V0Fitter::theLambdas
private

Definition at line 72 of file V0Fitter.h.

Referenced by fitAll(), and getLambdas().

double V0Fitter::tkChi2Cut
private

Definition at line 89 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::tkDCACut
private

Definition at line 99 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

int V0Fitter::tkNhitsCut
private

Definition at line 90 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

edm::EDGetTokenT<reco::BeamSpot> V0Fitter::token_beamSpot
private

Definition at line 105 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

edm::EDGetTokenT<reco::TrackCollection> V0Fitter::token_tracks
private

Definition at line 104 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

const TrackerGeometry* V0Fitter::trackerGeom
private

Definition at line 75 of file V0Fitter.h.

Referenced by fitAll().

bool V0Fitter::useRefTrax
private

Definition at line 79 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

edm::InputTag V0Fitter::vtxFitter
private

Definition at line 106 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::vtxSigCut
private

Definition at line 92 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().