CMS 3D CMS Logo

Public Member Functions | Private Attributes

InvariantMassAlgorithm Class Reference

#include <InvariantMassAlgorithm.h>

List of all members.

Public Member Functions

float getMinimumClusterDR (edm::Event &theEvent, const edm::EventSetup &theEventSetup, const reco::IsolatedTauTagInfoRef &tauRef, const math::XYZVector &cluster_3vec)
 InvariantMassAlgorithm ()
 InvariantMassAlgorithm (const edm::ParameterSet &parameters)
std::pair< double,
reco::TauMassTagInfo
tag (edm::Event &theEvent, const edm::EventSetup &theEventSetup, const reco::IsolatedTauTagInfoRef &tauRef, const edm::Handle< reco::BasicClusterCollection > &clus_handle)
 ~InvariantMassAlgorithm ()

Private Attributes

double cluster_jet_matching_cone
double cluster_track_matching_cone
double inv_mass_cut
double leading_trk_pt
double matching_cone
double signal_cone
TrackDetectorAssociatortrackAssociator_
TrackAssociatorParameters trackAssociatorParameters_

Detailed Description

Definition at line 26 of file InvariantMassAlgorithm.h.


Constructor & Destructor Documentation

InvariantMassAlgorithm::InvariantMassAlgorithm ( const edm::ParameterSet parameters)

Definition at line 23 of file InvariantMassAlgorithm.cc.

References edm::ParameterSet::getParameter().

{
  matching_cone               = parameters.getParameter<double>("MatchingCone");
  leading_trk_pt              = parameters.getParameter<double>("LeadingTrackPt");
  signal_cone                 = parameters.getParameter<double>("SignalCone");
  cluster_jet_matching_cone   = parameters.getParameter<double>("ClusterSelectionCone");
  cluster_track_matching_cone = parameters.getParameter<double>("ClusterTrackMatchingCone");
  inv_mass_cut                = parameters.getParameter<double>("InvariantMassCutoff");


// TrackAssociator parameters
   edm::ParameterSet tk_ass_pset = parameters.getParameter<edm::ParameterSet>("TrackAssociatorParameters");
   trackAssociatorParameters_.loadParameters( tk_ass_pset );

   trackAssociator_ = new TrackDetectorAssociator();
   trackAssociator_->useDefaultPropagator();
}
InvariantMassAlgorithm::InvariantMassAlgorithm ( )

Definition at line 11 of file InvariantMassAlgorithm.cc.

{
  matching_cone               = 0.1; // check with simone
  leading_trk_pt              = 1.0; // check with simone
  signal_cone                 = 0.07;// check with simone
  cluster_jet_matching_cone   = 0.4;
  cluster_track_matching_cone = 0.08;
  inv_mass_cut                = 2.5;
}
InvariantMassAlgorithm::~InvariantMassAlgorithm ( )

Definition at line 43 of file InvariantMassAlgorithm.cc.


Member Function Documentation

float InvariantMassAlgorithm::getMinimumClusterDR ( edm::Event theEvent,
const edm::EventSetup theEventSetup,
const reco::IsolatedTauTagInfoRef &  tauRef,
const math::XYZVector cluster_3vec 
)

Definition at line 86 of file InvariantMassAlgorithm.cc.

References edm::RefVector< C, T, F >::begin(), Geom::deltaR2(), edm::RefVector< C, T, F >::end(), info, metsig::jet, reco::LeafCandidate::px(), reco::LeafCandidate::py(), reco::LeafCandidate::pz(), testEve_cfg::tracks, and TrackDetMatchInfo::trkGlobPosAtEcal.

{

  trackAssociatorParameters_.useEcal = true;
  trackAssociatorParameters_.useHcal = false;
  trackAssociatorParameters_.useHO   = false;
  trackAssociatorParameters_.useMuon = false;
  trackAssociatorParameters_.dREcal  = 0.03;

  const TrackRefVector tracks = tauRef->allTracks();
  const Jet & jet = *(tauRef->jet()); 
  math::XYZVector jet3Vec(jet.px(),jet.py(),jet.pz());   
  float min_dR = 999.9; 
  float deltaR1 = ROOT::Math::VectorUtil::DeltaR(cluster_3vec, jet3Vec);
  if (deltaR1 > cluster_jet_matching_cone) return -1.0;

  for (edm::RefVector<reco::TrackCollection>::const_iterator it  = tracks.begin();
                                                             it != tracks.end(); it++) 
  {    
    TrackDetMatchInfo info;
    info = trackAssociator_->associate(theEvent, theEventSetup,
                                     trackAssociator_->getFreeTrajectoryState(theEventSetup, (*(*it))), 
                                     trackAssociatorParameters_);
    math::XYZVector track3Vec(info.trkGlobPosAtEcal.x(),
                              info.trkGlobPosAtEcal.y(),
                              info.trkGlobPosAtEcal.z());

    float deltaR2 = ROOT::Math::VectorUtil::DeltaR(cluster_3vec, track3Vec);
    if (deltaR2 < min_dR) min_dR = deltaR2;        
  }
  return min_dR;
}
pair< double, reco::TauMassTagInfo > InvariantMassAlgorithm::tag ( edm::Event theEvent,
const edm::EventSetup theEventSetup,
const reco::IsolatedTauTagInfoRef &  tauRef,
const edm::Handle< reco::BasicClusterCollection > &  clus_handle 
)

Definition at line 49 of file InvariantMassAlgorithm.cc.

References gather_cfg::cout, reco::TauMassTagInfo::discriminator(), reco::TauMassTagInfo::setIsolatedTauTag(), funct::sin(), reco::TauMassTagInfo::storeClusterTrackCollection(), x, detailsBasic3DVector::y, and z.

{

  TauMassTagInfo resultExtended;
  resultExtended.setIsolatedTauTag(tauRef);

  double discriminator = tauRef->discriminator();
  if (discriminator > 0.0) {

    int nSel = 0;
    for (size_t ic = 0; ic < clus_handle->size(); ic++) {
      math::XYZVector cluster3Vec((*clus_handle)[ic].x(),(*clus_handle)[ic].y(),(*clus_handle)[ic].z());
      BasicClusterRef cluster(clus_handle, ic);
      float et = (*clus_handle)[ic].energy() * sin(cluster3Vec.theta());       
      if ( (et < 0.04) || (et > 300.0) ) continue;    
      float delR = getMinimumClusterDR(theEvent, theEventSetup, tauRef, cluster3Vec);
      if (delR != -1.0) {
        nSel++;
        resultExtended.storeClusterTrackCollection(cluster, delR);
      }
    }
    if (0) cout << " Total Number of Clusters " << clus_handle->size() << " " << nSel << endl;
    
    discriminator = resultExtended.discriminator(matching_cone,
                                                 leading_trk_pt,
                                                 signal_cone,
                                                 cluster_track_matching_cone,
                                                 inv_mass_cut); 
  }
  return make_pair(discriminator, resultExtended); 
}

Member Data Documentation

Definition at line 50 of file InvariantMassAlgorithm.h.

Definition at line 51 of file InvariantMassAlgorithm.h.

Definition at line 52 of file InvariantMassAlgorithm.h.

Definition at line 48 of file InvariantMassAlgorithm.h.

Definition at line 47 of file InvariantMassAlgorithm.h.

Definition at line 49 of file InvariantMassAlgorithm.h.

Definition at line 54 of file InvariantMassAlgorithm.h.

Definition at line 55 of file InvariantMassAlgorithm.h.