CMS 3D CMS Logo

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