32 using namespace trigger;
37 srcTag_ (iConfig.getParameter<edm::
InputTag>(
"Src")),
39 previousCandTag_(iConfig.getParameter<edm::
InputTag>(
"PreviousCandTag")),
41 maxEta_ (iConfig.getParameter<double>(
"MaxEta")),
42 minPt_ (iConfig.getParameter<double>(
"MinPt")),
43 minPtPair_ (iConfig.getParameter<double>(
"MinPtPair")),
44 minInvMass_ (iConfig.getParameter<double>(
"MinInvMass")),
45 maxInvMass_ (iConfig.getParameter<double>(
"MaxInvMass")),
46 chargeOpt_ (iConfig.getParameter<int>(
"ChargeOpt"))
48 produces<VertexCollection>();
62 desc.
add<
double>(
"MaxEta",2.5);
63 desc.
add<
double>(
"MinPt",0.0);
64 desc.
add<
double>(
"MinPtPair",0.0);
65 desc.
add<
double>(
"MinInvMass",1.0);
66 desc.
add<
double>(
"MaxInvMass",20.0);
67 desc.
add<
int>(
"ChargeOpt",-1);
68 descriptions.
add(
"hltDisplacedmumuVtxProducer", desc);
86 double const MuMass = 0.106;
87 double const MuMass2 = MuMass*MuMass;
104 RecoChargedCandidateCollection::const_iterator cand1;
105 RecoChargedCandidateCollection::const_iterator cand2;
111 vector<RecoChargedCandidateRef> vPrevCands;
114 for (cand1=mucands->begin(); cand1!=mucands->end(); cand1++) {
116 LogDebug(
"HLTDisplacedMumuFilter") <<
" 1st muon in loop: q*pt= " << cand1->charge()*cand1->pt() <<
", eta= " << cand1->eta() <<
", hits= " << tk1->numberOfValidHits();
122 if (fabs(cand1->eta())>
maxEta_)
continue;
123 if (cand1->pt() <
minPt_)
continue;
125 cand2 = cand1; cand2++;
126 for (; cand2!=mucands->end(); cand2++) {
130 LogDebug(
"HLTDisplacedmumuVtxProducer") <<
" 2nd muon in loop: q*pt= " << cand2->charge()*cand2->pt() <<
", eta= " << cand2->eta() <<
", hits= " << tk2->numberOfValidHits() <<
", d0= " << tk2->d0();
135 if (fabs(cand2->eta())>
maxEta_)
continue;
136 if (cand2->pt() <
minPt_)
continue;
140 if (cand1->charge()*cand2->charge()>0)
continue;
142 if (cand1->charge()*cand2->charge()<0)
continue;
146 e1 =
sqrt(cand1->momentum().Mag2()+MuMass2);
147 e2 =
sqrt(cand2->momentum().Mag2()+MuMass2);
155 double invmass =
abs(p.mass());
156 LogDebug(
"HLTDisplacedMumuFilter") <<
" ... 1-2 invmass= " << invmass;
162 vector<TransientTrack> t_tks;
165 t_tks.push_back(ttkp1);
166 t_tks.push_back(ttkp2);
169 if (t_tks.size()!=2)
continue;
179 vertexCollection->push_back(vertex);
182 iEvent.
put(vertexCollection);
189 for (
unsigned int i=0;
i<refVect.size();
i++) {
190 if ( refVect[
i]->get<TrackRef>() == trackref ) {
bool checkPreviousCand(const reco::TrackRef &trackref, std::vector< reco::RecoChargedCandidateRef > &ref2)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
HLTDisplacedmumuVtxProducer(const edm::ParameterSet &)
virtual CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const
edm::EDGetTokenT< reco::RecoChargedCandidateCollection > srcToken_
std::vector< Vertex > VertexCollection
collection of Vertex objects
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > previousCandToken_
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
~HLTDisplacedmumuVtxProducer()
Abs< T >::type abs(const T &t)
virtual void produce(edm::Event &, const edm::EventSetup &)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
std::vector< RecoChargedCandidate > RecoChargedCandidateCollection
collectin of RecoChargedCandidate objects
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
T const * get() const
Returns C++ pointer to the item.
math::XYZTLorentzVector LorentzVector
Lorentz vector.
math::PtEtaPhiELorentzVectorF LorentzVector