CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTmumutkVtxProducer.cc
Go to the documentation of this file.
1 #include <algorithm>
2 #include <cmath>
3 #include <vector>
4 
8 
21 
24 
25 using namespace edm;
26 using namespace reco;
27 using namespace std;
28 using namespace trigger;
29 
30 // ----------------------------------------------------------------------
32  muCandTag_ (iConfig.getParameter<edm::InputTag>("MuCand")),
33  muCandToken_(consumes<reco::RecoChargedCandidateCollection>(muCandTag_)),
34  trkCandTag_ (iConfig.getParameter<edm::InputTag>("TrackCand")),
35  trkCandToken_(consumes<reco::RecoChargedCandidateCollection>(trkCandTag_)),
36  previousCandTag_(iConfig.getParameter<edm::InputTag>("PreviousCandTag")),
37  previousCandToken_(consumes<trigger::TriggerFilterObjectWithRefs>(previousCandTag_)),
38  mfName_(iConfig.getParameter<std::string>("SimpleMagneticField")),
39  thirdTrackMass_(iConfig.getParameter<double>("ThirdTrackMass")),
40  maxEta_(iConfig.getParameter<double>("MaxEta")),
41  minPt_(iConfig.getParameter<double>("MinPt")),
42  minInvMass_(iConfig.getParameter<double>("MinInvMass")),
43  maxInvMass_(iConfig.getParameter<double>("MaxInvMass")),
44  minD0Significance_(iConfig.getParameter<double>("MinD0Significance")),
45  beamSpotTag_ (iConfig.getParameter<edm::InputTag> ("BeamSpotTag")),
46  beamSpotToken_(consumes<reco::BeamSpot>(beamSpotTag_))
47 {
48  produces<VertexCollection>();
49 }
50 
51 // ----------------------------------------------------------------------
53 
54 void
57  desc.add<edm::InputTag>("MuCand",edm::InputTag("hltMuTracks"));
58  desc.add<edm::InputTag>("TrackCand",edm::InputTag("hltMumukAllConeTracks"));
59  desc.add<edm::InputTag>("PreviousCandTag",edm::InputTag("hltDisplacedmumuFilterDoubleMu4Jpsi"));
60  desc.add<std::string>("SimpleMagneticField","");
61  desc.add<double>("ThirdTrackMass",0.493677);
62  desc.add<double>("MaxEta",2.5);
63  desc.add<double>("MinPt",3.0);
64  desc.add<double>("MinInvMass",0.0);
65  desc.add<double>("MaxInvMass",99999.);
66  desc.add<double>("MinD0Significance",0.0);
67  desc.add<edm::InputTag>("BeamSpotTag",edm::InputTag("hltOfflineBeamSpot"));
68  descriptions.add("HLTmumutkVtxProducer",desc);
69 }
70 
71 // ----------------------------------------------------------------------
73 {
74  const double MuMass(0.106);
75  const double MuMass2(MuMass*MuMass);
76  const double thirdTrackMass2(thirdTrackMass_*thirdTrackMass_);
77 
78  // get hold of muon trks
80  iEvent.getByToken(muCandToken_,mucands);
81 
82  //get the transient track builder:
84  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theB);
85 
86  //get the beamspot position
87  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
88  iEvent.getByToken(beamSpotToken_,recoBeamSpotHandle);
89 
90  //get the b field
91  ESHandle<MagneticField> bFieldHandle;
92  iSetup.get<IdealMagneticFieldRecord>().get(mfName_, bFieldHandle);
93  const MagneticField* magField = bFieldHandle.product();
94  TSCBLBuilderNoMaterial blsBuilder;
95 
96  // get track candidates around displaced muons
98  iEvent.getByToken(trkCandToken_,trkcands);
99 
100  auto_ptr<VertexCollection> vertexCollection(new VertexCollection());
101 
102  // Ref to Candidate object to be recorded in filter object
105  RecoChargedCandidateRef refTrk;
106 
107  double e1,e2,e3;
109 
110  if ( mucands->size() < 2 ) return;
111  if ( trkcands->size() < 1 ) return;
112 
113  RecoChargedCandidateCollection::const_iterator mucand1;
114  RecoChargedCandidateCollection::const_iterator mucand2;
115  RecoChargedCandidateCollection::const_iterator trkcand;
116 
117  // get the objects passing the previous filter
119  iEvent.getByToken(previousCandToken_,previousCands);
120 
121  vector<RecoChargedCandidateRef> vPrevCands;
122  previousCands->getObjects(TriggerMuon,vPrevCands);
123 
124  for (mucand1=mucands->begin(); mucand1!=mucands->end(); ++mucand1) {
125  TrackRef trk1 = mucand1->get<TrackRef>();
126  LogDebug("HLTmumutkVtxProducer") << " 1st muon: q*pt= " << trk1->charge()*trk1->pt()
127  << ", eta= " << trk1->eta()
128  << ", hits= " << trk1->numberOfValidHits();
129 
130  //first check if this muon passed the previous filter
131  if( ! checkPreviousCand( trk1, vPrevCands) ) continue;
132 
133  // eta and pt cut
134  if (fabs(trk1->eta()) > maxEta_) continue;
135  if (trk1->pt() < minPt_ ) continue;
136 
137  mucand2 = mucand1; ++mucand2;
138  for (; mucand2!=mucands->end(); mucand2++) {
139  TrackRef trk2 = mucand2->get<TrackRef>();
140 
141  LogDebug("HLTDisplacedMumukFilter") << " 2nd muon: q*pt= " << trk2->charge()*trk2->pt()
142  << ", eta= " << trk2->eta()
143  << ", hits= " << trk2->numberOfValidHits();
144 
145  //first check if this muon passed the previous filter
146  if( ! checkPreviousCand( trk2, vPrevCands) ) continue;
147  // eta and pt cut
148  if (fabs(trk2->eta()) > maxEta_) continue;
149  if (trk2->pt() < minPt_ ) continue;
150 
151  //loop on track collection
152  for ( trkcand = trkcands->begin(); trkcand !=trkcands->end(); ++trkcand) {
153  TrackRef trk3 = trkcand->get<TrackRef>();
154  if( overlap( trk1, trk3) ) continue;
155  if( overlap( trk2, trk3) ) continue;
156 
157  LogDebug("HLTDisplacedMumukFilter") << " 3rd track: q*pt= " << trk3->charge()*trk3->pt()
158  << ", eta= " << trk3->eta()
159  << ", hits= " << trk3->numberOfValidHits();
160 
161  // eta and pt cut
162  if (fabs(trk3->eta()) > maxEta_) continue;
163  if (trk3->pt() < minPt_ ) continue;
164 
165  // Combined system
166  e1 = sqrt(trk1->momentum().Mag2() + MuMass2 );
167  e2 = sqrt(trk2->momentum().Mag2() + MuMass2 );
168  e3 = sqrt(trk3->momentum().Mag2() + thirdTrackMass2);
169 
170  p1 = Particle::LorentzVector(trk1->px(),trk1->py(),trk1->pz(),e1);
171  p2 = Particle::LorentzVector(trk2->px(),trk2->py(),trk2->pz(),e2);
172  p3 = Particle::LorentzVector(trk3->px(),trk3->py(),trk3->pz(),e3);
173 
174  p = p1+p2+p3;
175 
176  //invariant mass cut
177  double invmass = abs(p.mass());
178  LogDebug("HLTDisplacedMumukFilter") << " Invmass= " << invmass;
179  if (invmass<minInvMass_) continue;
180  if (invmass>maxInvMass_) continue;
181 
182  // do the vertex fit
183  vector<TransientTrack> t_tks;
184  t_tks.push_back((*theB).build(&trk1));
185  t_tks.push_back((*theB).build(&trk2));
186  t_tks.push_back((*theB).build(&trk3));
187  if (t_tks.size()!=3) continue;
188 
189  FreeTrajectoryState InitialFTS = initialFreeState(*trk3, magField);
190  TrajectoryStateClosestToBeamLine tscb( blsBuilder(InitialFTS, *recoBeamSpotHandle) );
191  double d0sig = tscb.transverseImpactParameter().significance();
192  if (d0sig < minD0Significance_) continue;
193 
194  KalmanVertexFitter kvf;
195  TransientVertex tv = kvf.vertex(t_tks);
196  if (!tv.isValid()) continue;
197  Vertex vertex = tv;
198 
199  // put vertex in the event
200  vertexCollection->push_back(vertex);
201  }
202  }
203  }
204  iEvent.put(vertexCollection);
205 }
206 
208 {
209  Basic3DVector<float> pos( tk.vertex());
210  GlobalPoint gpos( pos);
211  Basic3DVector<float> mom( tk.momentum());
212  GlobalVector gmom( mom);
213  GlobalTrajectoryParameters par( gpos, gmom, tk.charge(), field);
215  return FreeTrajectoryState( par, err);
216 }
217 
218 bool HLTmumutkVtxProducer::overlap(const TrackRef& trackref1, const TrackRef& trackref2){
219  double eps(1.44e-4);
220  if (deltaR(trackref1->eta(), trackref1->phi(),trackref2->eta(), trackref2->phi()) < eps) return 1;
221  return 0;
222 }
223 
224 bool HLTmumutkVtxProducer::checkPreviousCand(const TrackRef& trackref, vector<RecoChargedCandidateRef> & refVect){
225  bool ok=false;
226  for (unsigned int i=0; i<refVect.size(); i++) {
227  if ( refVect[i]->get<TrackRef>() == trackref ) {
228  ok=true;
229  break;
230  }
231  }
232  return ok;
233 }
#define LogDebug(id)
edm::EDGetTokenT< reco::RecoChargedCandidateCollection > muCandToken_
static FreeTrajectoryState initialFreeState(const reco::Track &, const MagneticField *)
const std::string mfName_
int i
Definition: DBlmapReader.cc:9
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > previousCandToken_
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:148
edm::EDGetTokenT< reco::RecoChargedCandidateCollection > trkCandToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
virtual CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::EDGetTokenT< reco::BeamSpot > beamSpotToken_
tuple vertexCollection
CovarianceMatrix covariance() const
return track covariance matrix
Definition: TrackBase.h:180
bool overlap(const reco::TrackRef &trackref1, const reco::TrackRef &trackref2)
int iEvent
Definition: GenABIO.cc:230
HLTmumutkVtxProducer(const edm::ParameterSet &)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
T sqrt(T t)
Definition: SSEVec.h:48
bool checkPreviousCand(const reco::TrackRef &trackref, std::vector< reco::RecoChargedCandidateRef > &ref2)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ParameterDescriptionBase * add(U const &iLabel, T const &value)
double p2[4]
Definition: TauolaWrapper.h:90
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
const Point & vertex() const
reference point on the track. This method is DEPRECATED, please use referencePoint() instead ...
Definition: TrackBase.h:154
std::vector< RecoChargedCandidate > RecoChargedCandidateCollection
collectin of RecoChargedCandidate objects
Float e1
Definition: deltaR.h:22
double significance() const
Definition: Measurement1D.h:32
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Float e2
Definition: deltaR.h:23
double p1[4]
Definition: TauolaWrapper.h:89
int charge() const
track electric charge
Definition: TrackBase.h:111
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:242
virtual void produce(edm::Event &, const edm::EventSetup &)
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:27
bool isValid() const
math::PtEtaPhiELorentzVectorF LorentzVector
double p3[4]
Definition: TauolaWrapper.h:91
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)