CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

HLTDisplacedmumuVtxProducer Class Reference

#include <HLTDisplacedmumuVtxProducer.h>

Inheritance diagram for HLTDisplacedmumuVtxProducer:
edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

virtual void beginJob ()
virtual void endJob ()
 HLTDisplacedmumuVtxProducer (const edm::ParameterSet &)
virtual void produce (edm::Event &, const edm::EventSetup &)
 ~HLTDisplacedmumuVtxProducer ()

Private Member Functions

bool checkPreviousCand (const reco::TrackRef &trackref, std::vector< reco::RecoChargedCandidateRef > &ref2)

Private Attributes

int chargeOpt_
double maxEta_
double maxInvMass_
double minInvMass_
double minPt_
double minPtPair_
edm::InputTag previousCandTag_
edm::InputTag src_

Detailed Description

Definition at line 28 of file HLTDisplacedmumuVtxProducer.h.


Constructor & Destructor Documentation

HLTDisplacedmumuVtxProducer::HLTDisplacedmumuVtxProducer ( const edm::ParameterSet iConfig) [explicit]

Definition at line 35 of file HLTDisplacedmumuVtxProducer.cc.

                                                                                      : 
        src_ (iConfig.getParameter<edm::InputTag>("Src")),
        previousCandTag_(iConfig.getParameter<edm::InputTag>("PreviousCandTag")),
        maxEta_ (iConfig.getParameter<double>("MaxEta")),
        minPt_ (iConfig.getParameter<double>("MinPt")),
        minPtPair_ (iConfig.getParameter<double>("MinPtPair")),
        minInvMass_ (iConfig.getParameter<double>("MinInvMass")),
        maxInvMass_ (iConfig.getParameter<double>("MaxInvMass")),
        chargeOpt_ (iConfig.getParameter<int>("ChargeOpt"))
{
        produces<VertexCollection>();
}
HLTDisplacedmumuVtxProducer::~HLTDisplacedmumuVtxProducer ( )

Definition at line 49 of file HLTDisplacedmumuVtxProducer.cc.

{

}

Member Function Documentation

void HLTDisplacedmumuVtxProducer::beginJob ( void  ) [virtual]

Reimplemented from edm::EDProducer.

Definition at line 56 of file HLTDisplacedmumuVtxProducer.cc.

{

}
bool HLTDisplacedmumuVtxProducer::checkPreviousCand ( const reco::TrackRef trackref,
std::vector< reco::RecoChargedCandidateRef > &  ref2 
) [private]

Definition at line 188 of file HLTDisplacedmumumuVtxProducer.cc.

References i, and convertSQLiteXML::ok.

Referenced by produce().

                                                                                                                        {
  bool ok=false;
  for (unsigned int i=0; i<refVect.size(); i++) {
    if ( refVect[i]->get<TrackRef>() == trackref ) {
      ok=true;
      break;
    }
  }
  return ok;
}
void HLTDisplacedmumuVtxProducer::endJob ( void  ) [virtual]

Reimplemented from edm::EDProducer.

Definition at line 62 of file HLTDisplacedmumuVtxProducer.cc.

{
        
}
void HLTDisplacedmumuVtxProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
) [virtual]

Implements edm::EDProducer.

Definition at line 68 of file HLTDisplacedmumuVtxProducer.cc.

References abs, chargeOpt_, checkPreviousCand(), edm::EventSetup::get(), edm::Ref< C, T, F >::get(), edm::Event::getByLabel(), TransientVertex::isValid(), LogDebug, maxEta_, maxInvMass_, minInvMass_, minPt_, minPtPair_, AlCaHLTBitMon_ParallelJobs::p, p1, p2, previousCandTag_, edm::Event::put(), mathSSE::sqrt(), src_, trigger::TriggerMuon, KalmanVertexFitter::vertex(), and GoodVertex_cfg::vertexCollection.

{
        double const MuMass = 0.106;
        double const MuMass2 = MuMass*MuMass;
        
        
        // get hold of muon trks
        Handle<RecoChargedCandidateCollection> mucands;
        iEvent.getByLabel (src_,mucands);
        
        //get the transient track builder:
        edm::ESHandle<TransientTrackBuilder> theB;
        iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theB);

        std::auto_ptr<VertexCollection> vertexCollection(new VertexCollection());

        // look at all mucands,  check cuts and make vertices
        double e1,e2;
        Particle::LorentzVector p,p1,p2;
        
        RecoChargedCandidateCollection::const_iterator cand1;
        RecoChargedCandidateCollection::const_iterator cand2;
        
        // get the objects passing the previous filter
        Handle<TriggerFilterObjectWithRefs> previousCands;
        iEvent.getByLabel (previousCandTag_,previousCands);

        vector<RecoChargedCandidateRef> vPrevCands;
        previousCands->getObjects(TriggerMuon,vPrevCands);

        for (cand1=mucands->begin(); cand1!=mucands->end(); cand1++) {
               TrackRef tk1 = cand1->get<TrackRef>();
               LogDebug("HLTDisplacedMumuFilter") << " 1st muon in loop: q*pt= " << cand1->charge()*cand1->pt() << ", eta= " << cand1->eta() << ", hits= " << tk1->numberOfValidHits();
             
               //first check if this muon passed the previous filter
               if( ! checkPreviousCand( tk1, vPrevCands) ) continue;
        
               // cuts
               if (fabs(cand1->eta())>maxEta_) continue;
               if (cand1->pt() < minPt_) continue;
              
               cand2 = cand1; cand2++;
               for (; cand2!=mucands->end(); cand2++) {
                         TrackRef tk2 = cand2->get<TrackRef>();
                 
                         // eta cut
                         LogDebug("HLTDisplacedmumuVtxProducer") << " 2nd muon in loop: q*pt= " << cand2->charge()*cand2->pt() << ", eta= " << cand2->eta() << ", hits= " << tk2->numberOfValidHits() << ", d0= " << tk2->d0();
                         //first check if this muon passed the previous filter
                         if( ! checkPreviousCand( tk2, vPrevCands) ) continue;
                         
                         // cuts
                         if (fabs(cand2->eta())>maxEta_) continue;
                         if (cand2->pt() < minPt_) continue;
                         
                         // opposite sign or same sign
                         if (chargeOpt_<0) {
                           if (cand1->charge()*cand2->charge()>0) continue;
                         } else if (chargeOpt_>0) {
                           if (cand1->charge()*cand2->charge()<0) continue;
                         }
                         
                         // Combined dimuon system
                         e1 = sqrt(cand1->momentum().Mag2()+MuMass2);
                         e2 = sqrt(cand2->momentum().Mag2()+MuMass2);
                         p1 = Particle::LorentzVector(cand1->px(),cand1->py(),cand1->pz(),e1);
                         p2 = Particle::LorentzVector(cand2->px(),cand2->py(),cand2->pz(),e2);
                         p = p1+p2;
                         
                         
                         if (p.pt()<minPtPair_) continue;
                         
                         double invmass = abs(p.mass());
                         LogDebug("HLTDisplacedMumuFilter") << " ... 1-2 invmass= " << invmass;
                         
                         if (invmass<minInvMass_) continue;
                         if (invmass>maxInvMass_) continue;
                         
                         // do the vertex fit
                         vector<TransientTrack> t_tks;
                         TransientTrack ttkp1 = (*theB).build(&tk1);
                         TransientTrack ttkp2 = (*theB).build(&tk2);
                         t_tks.push_back(ttkp1);
                         t_tks.push_back(ttkp2);
                         
                
                         if (t_tks.size()!=2) continue;
                        
                         KalmanVertexFitter kvf;
                         TransientVertex tv = kvf.vertex(t_tks);

                         if (!tv.isValid()) continue;
                         
                         Vertex vertex = tv;

                         // put vertex in the event
                         vertexCollection->push_back(vertex);
               }
        }
        iEvent.put(vertexCollection);
}

Member Data Documentation

Definition at line 47 of file HLTDisplacedmumuVtxProducer.h.

Referenced by produce().

Definition at line 42 of file HLTDisplacedmumuVtxProducer.h.

Referenced by produce().

Definition at line 46 of file HLTDisplacedmumuVtxProducer.h.

Referenced by produce().

Definition at line 45 of file HLTDisplacedmumuVtxProducer.h.

Referenced by produce().

Definition at line 43 of file HLTDisplacedmumuVtxProducer.h.

Referenced by produce().

Definition at line 44 of file HLTDisplacedmumuVtxProducer.h.

Referenced by produce().

Definition at line 41 of file HLTDisplacedmumuVtxProducer.h.

Referenced by produce().

Definition at line 40 of file HLTDisplacedmumuVtxProducer.h.

Referenced by produce().