CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes

FML1PtSmearer Class Reference

#include <FML1PtSmearer.h>

List of all members.

Public Member Functions

 FML1PtSmearer (const RandomEngine *engine)
 Constructor (read in the parametrizations from a data file)
bool smear (SimpleL1MuGMTCand *)
 smear the transverse momentum of a SimplL1MuGMTCand
 ~FML1PtSmearer ()
 Destructor.

Private Member Functions

float ChargeMisIdent (int ieta, double pt)
int IndexOfPtgen (float pt)

Private Attributes

const RandomEnginerandom
float resolution [DIMRES]

Static Private Attributes

static const int DIMRES = 3*NPT*NPTL1
static const int NPT = 136
static const int NPTL1 = 31

Detailed Description

Class to deal with the 'smearing' of the L1 muon transverse momentum. The output momentum is generated according to the probablility that a MC muon with the same pt leads to that (discrete) momentum value in the GMT.

Author:
Andrea Perrotta Date: 26/11/2004 05/09/2006

Definition at line 16 of file FML1PtSmearer.h.


Constructor & Destructor Documentation

FML1PtSmearer::FML1PtSmearer ( const RandomEngine engine)

Constructor (read in the parametrizations from a data file)

Definition at line 16 of file FML1PtSmearer.cc.

References gather_cfg::cout, connectstrParser::f1, alignmentValidation::fname, i, j, NPT, NPTL1, path(), and resolution.

  : random(engine) {

  string fname = "FastSimulation/Muons/data/resolutionL1.data";
  //   std::string path(std::getenv("CMSSW_SEARCH__PATH"));
  std::string path(getenv("CMSSW_SEARCH_PATH"));
  FileInPath f1(path,fname);
  if ( f1() == 0) {
    std::cout << "File " << fname << " not found in " << path << std::endl;
    throw Genexception(" resolution list not found for FastMuonLvl1Trigger.");
  } else {
    // The following should be on LogInfo
    //cout << "Reading " << f1.name() << std::endl;
  }
  std::ifstream & listfile = *f1();
   
  int ind;
  float num=0.;
  for (int ieta=0; ieta<3; ieta++) {
    for (int i=0; i<NPT; i++) {
      float sum = 0.;
      ind = i*NPTL1 + ieta*(NPT*NPTL1);
      for (int j=0; j<NPTL1; j++) {
        listfile >> num;
        sum += num;
        resolution[ind]=sum; 
        ind++;
      }
      ind = i*NPTL1 + ieta*(NPT*NPTL1);
      for (int j=0; j<NPTL1 ; j++){
        resolution[ind] /= sum;
        ind++;
      }
    }
  }
  
}
FML1PtSmearer::~FML1PtSmearer ( )

Destructor.

Definition at line 54 of file FML1PtSmearer.cc.

{}

Member Function Documentation

float FML1PtSmearer::ChargeMisIdent ( int  ieta,
double  pt 
) [inline, private]

Definition at line 40 of file FML1PtSmearer.h.

References ExpressReco_HICollisions_FallBack::pt.

Referenced by smear().

                                                    {
    float df=0.;
    switch (ieta) {
    case 0: 
      df =  2.16909e-03  + 1.95708e-04*pt ;
      break;
    case 1:  
      if (pt>500.) pt = 500.;
      df = 1.00445e-02 + 1.15253e-03*pt - 7.73819e-07*pt*pt ;
      break;
    case 2:
      if (pt>200.) pt = 200.;
      df = 5.00580e-02 + 5.88949e-03*pt - 2.98100e-05*pt*pt + 5.02454e-08*pt*pt*pt ;
      break;
    }
    return (df<0.5? df: 0.5) ;
  }
int FML1PtSmearer::IndexOfPtgen ( float  pt) [private]

Definition at line 100 of file FML1PtSmearer.cc.

References i, and NPT.

Referenced by smear().

                                        {

  static float vecpt[NPT] = {
                1.,   2.,   3.,   4.,   5.,   6.,   7.,   8.,   9.,
         10.,  11.,  12.,  13.,  14.,  15.,  16.,  17.,  18.,  19.,
         20.,  21.,  22.,  23.,  24.,  25.,  26.,  27.,  28.,  29.,
         30.,  31.,  32.,  33.,  34.,  35.,  36.,  37.,  38.,  39.,
         40.,  41.,  42.,  43.,  44.,  45.,  46.,  47.,  48.,  49.,
         50.,  51.,  52.,  53.,  54.,  55.,  56.,  57.,  58.,  59.,
         60.,  61.,  62.,  63.,  64.,  65.,  66.,  67.,  68.,  69.,
         70.,  71.,  72.,  73.,  74.,  75.,  76.,  77.,  78.,  79.,
         80.,  81.,  82.,  83.,  84.,  85.,  86.,  87.,  88.,  89.,
         90.,  91.,  92.,  93.,  94.,  95.,  96.,  97.,  98.,  99.,
        100., 110., 120., 130., 140., 150., 160., 170., 180., 190.,
        200., 210., 220., 230., 240., 250., 260., 270., 280., 290.,
        300., 310., 320., 330., 340., 350., 360., 370., 380., 390.,
        400., 500., 600., 700., 800., 900., 1000. };

  for (int i=0; i<NPT; i++) {
    if (pt<vecpt[i]) return i;
  }
  return NPT;
}
bool FML1PtSmearer::smear ( SimpleL1MuGMTCand aMuon)

smear the transverse momentum of a SimplL1MuGMTCand

Definition at line 58 of file FML1PtSmearer.cc.

References SimpleL1MuGMTCand::charge(), ChargeMisIdent(), cmsDriverOptions::counter, RandomEngine::flatShoot(), SimpleL1MuGMTCand::getMomentum(), IndexOfPtgen(), NPT, NPTL1, reco::tau::disc::Pt(), SimpleL1MuGMTCand::ptScale, random, resolution, SimpleL1MuGMTCand::setCharge(), L1MuGMTCand::setPtPacked(), L1MuGMTCand::setPtValue(), and SimpleL1MuGMTCand::smearedPt().

Referenced by ParamL3MuonProducer::produce(), and FastL1MuonProducer::produce().

                                                   {
  // get and smear the pt
  double Pt=aMuon->smearedPt();
  double AbsEta=fabs(aMuon->getMomentum().eta());

  bool StatusCode=true;  
  
  if (AbsEta>2.40 || Pt<3) { StatusCode=false;}

  if (StatusCode){

    int ieta = 0;
    if      (AbsEta<1.04) ieta = 0 ;
    else if (AbsEta<2.07) ieta = 1 ;
    else if (AbsEta<2.40) ieta = 2 ;

    int counter = 0;
    int ipt = IndexOfPtgen(Pt);
    int ind = counter + ipt*NPTL1 + ieta*(NPT*NPTL1);
    double prob = random->flatShoot();
    while ( (prob > resolution[ind]) && counter<NPTL1 ){ 
      counter++;
      ind++;
    }
    counter++;
    aMuon->setPtPacked(counter & 31);
    aMuon->setPtValue(SimpleL1MuGMTCand::ptScale[counter]);

    prob = random->flatShoot();
    if (prob <= ChargeMisIdent(ieta,Pt)) {
      int OldCharge = aMuon->charge(); 
      int NewCharge = -OldCharge;
      aMuon->setCharge(NewCharge);
    } 

  }

  return StatusCode;

}  

Member Data Documentation

const int FML1PtSmearer::DIMRES = 3*NPT*NPTL1 [static, private]

Definition at line 37 of file FML1PtSmearer.h.

const int FML1PtSmearer::NPT = 136 [static, private]

Definition at line 36 of file FML1PtSmearer.h.

Referenced by FML1PtSmearer(), IndexOfPtgen(), and smear().

const int FML1PtSmearer::NPTL1 = 31 [static, private]

Definition at line 35 of file FML1PtSmearer.h.

Referenced by FML1PtSmearer(), and smear().

Definition at line 31 of file FML1PtSmearer.h.

Referenced by smear().

Definition at line 38 of file FML1PtSmearer.h.

Referenced by FML1PtSmearer(), and smear().