CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
FML1PtSmearer Class Reference

#include <FML1PtSmearer.h>

Public Member Functions

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

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.

17  : random(engine) {
18 
19  string fname = "FastSimulation/Muons/data/resolutionL1.data";
20  // std::string path(std::getenv("CMSSW_SEARCH__PATH"));
21  std::string path(getenv("CMSSW_SEARCH_PATH"));
22  FileInPath f1(path,fname);
23  if ( f1() == 0) {
24  std::cout << "File " << fname << " not found in " << path << std::endl;
25  throw Genexception(" resolution list not found for FastMuonLvl1Trigger.");
26  } else {
27  // The following should be on LogInfo
28  //cout << "Reading " << f1.name() << std::endl;
29  }
30  std::ifstream & listfile = *f1();
31 
32  int ind;
33  float num=0.;
34  for (int ieta=0; ieta<3; ieta++) {
35  for (int i=0; i<NPT; i++) {
36  float sum = 0.;
37  ind = i*NPTL1 + ieta*(NPT*NPTL1);
38  for (int j=0; j<NPTL1; j++) {
39  listfile >> num;
40  sum += num;
41  resolution[ind]=sum;
42  ind++;
43  }
44  ind = i*NPTL1 + ieta*(NPT*NPTL1);
45  for (int j=0; j<NPTL1 ; j++){
46  resolution[ind] /= sum;
47  ind++;
48  }
49  }
50  }
51 
52 }
int i
Definition: DBlmapReader.cc:9
static const int NPTL1
Definition: FML1PtSmearer.h:35
const RandomEngine * random
Definition: FML1PtSmearer.h:31
static const int NPT
Definition: FML1PtSmearer.h:36
float resolution[DIMRES]
Definition: FML1PtSmearer.h:38
int path() const
Definition: HLTadd.h:3
int j
Definition: DBlmapReader.cc:9
long long int num
Definition: procUtils.cc:71
string fname
main script
tuple cout
Definition: gather_cfg.py:41
FML1PtSmearer::~FML1PtSmearer ( )

Destructor.

Definition at line 54 of file FML1PtSmearer.cc.

54 {}

Member Function Documentation

float FML1PtSmearer::ChargeMisIdent ( int  ieta,
double  pt 
)
inlineprivate

Definition at line 40 of file FML1PtSmearer.h.

References ExpressReco_HICollisions_FallBack::pt.

Referenced by smear().

40  {
41  float df=0.;
42  switch (ieta) {
43  case 0:
44  df = 2.16909e-03 + 1.95708e-04*pt ;
45  break;
46  case 1:
47  if (pt>500.) pt = 500.;
48  df = 1.00445e-02 + 1.15253e-03*pt - 7.73819e-07*pt*pt ;
49  break;
50  case 2:
51  if (pt>200.) pt = 200.;
52  df = 5.00580e-02 + 5.88949e-03*pt - 2.98100e-05*pt*pt + 5.02454e-08*pt*pt*pt ;
53  break;
54  }
55  return (df<0.5? df: 0.5) ;
56  }
int FML1PtSmearer::IndexOfPtgen ( float  pt)
private

Definition at line 100 of file FML1PtSmearer.cc.

References i, and NPT.

Referenced by smear().

100  {
101 
102  static float vecpt[NPT] = {
103  1., 2., 3., 4., 5., 6., 7., 8., 9.,
104  10., 11., 12., 13., 14., 15., 16., 17., 18., 19.,
105  20., 21., 22., 23., 24., 25., 26., 27., 28., 29.,
106  30., 31., 32., 33., 34., 35., 36., 37., 38., 39.,
107  40., 41., 42., 43., 44., 45., 46., 47., 48., 49.,
108  50., 51., 52., 53., 54., 55., 56., 57., 58., 59.,
109  60., 61., 62., 63., 64., 65., 66., 67., 68., 69.,
110  70., 71., 72., 73., 74., 75., 76., 77., 78., 79.,
111  80., 81., 82., 83., 84., 85., 86., 87., 88., 89.,
112  90., 91., 92., 93., 94., 95., 96., 97., 98., 99.,
113  100., 110., 120., 130., 140., 150., 160., 170., 180., 190.,
114  200., 210., 220., 230., 240., 250., 260., 270., 280., 290.,
115  300., 310., 320., 330., 340., 350., 360., 370., 380., 390.,
116  400., 500., 600., 700., 800., 900., 1000. };
117 
118  for (int i=0; i<NPT; i++) {
119  if (pt<vecpt[i]) return i;
120  }
121  return NPT;
122 }
int i
Definition: DBlmapReader.cc:9
static const int NPT
Definition: FML1PtSmearer.h:36
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 FastL1MuonProducer::produce(), and ParamL3MuonProducer::produce().

58  {
59  // get and smear the pt
60  double Pt=aMuon->smearedPt();
61  double AbsEta=fabs(aMuon->getMomentum().eta());
62 
63  bool StatusCode=true;
64 
65  if (AbsEta>2.40 || Pt<3) { StatusCode=false;}
66 
67  if (StatusCode){
68 
69  int ieta = 0;
70  if (AbsEta<1.04) ieta = 0 ;
71  else if (AbsEta<2.07) ieta = 1 ;
72  else if (AbsEta<2.40) ieta = 2 ;
73 
74  int counter = 0;
75  int ipt = IndexOfPtgen(Pt);
76  int ind = counter + ipt*NPTL1 + ieta*(NPT*NPTL1);
77  double prob = random->flatShoot();
78  while ( (prob > resolution[ind]) && counter<NPTL1 ){
79  counter++;
80  ind++;
81  }
82  counter++;
83  aMuon->setPtPacked(counter & 31);
84  aMuon->setPtValue(SimpleL1MuGMTCand::ptScale[counter]);
85 
86  prob = random->flatShoot();
87  if (prob <= ChargeMisIdent(ieta,Pt)) {
88  int OldCharge = aMuon->charge();
89  int NewCharge = -OldCharge;
90  aMuon->setCharge(NewCharge);
91  }
92 
93  }
94 
95  return StatusCode;
96 
97 }
void setPtPacked(unsigned pt)
set packed pt-code of muon candidate
Definition: L1MuGMTCand.h:156
static const int NPTL1
Definition: FML1PtSmearer.h:35
const RandomEngine * random
Definition: FML1PtSmearer.h:31
static const int NPT
Definition: FML1PtSmearer.h:36
float resolution[DIMRES]
Definition: FML1PtSmearer.h:38
int charge() const
get charge
static const float ptScale[32]
void setCharge(int charge)
set charge and packed code of muon candidate
int IndexOfPtgen(float pt)
float ChargeMisIdent(int ieta, double pt)
Definition: FML1PtSmearer.h:40
float smearedPt() const
return the smeared L1 Pt value before discretization in 32-bit
const LorentzVector getMomentum() const
double flatShoot(double xmin=0.0, double xmax=1.0) const
Definition: RandomEngine.h:30
void setPtValue(float ptVal)
Set Pt Value.
Definition: L1MuGMTCand.h:182

Member Data Documentation

const int FML1PtSmearer::DIMRES =3*NPT*NPTL1
staticprivate

Definition at line 37 of file FML1PtSmearer.h.

const int FML1PtSmearer::NPT =136
staticprivate

Definition at line 36 of file FML1PtSmearer.h.

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

const int FML1PtSmearer::NPTL1 =31
staticprivate

Definition at line 35 of file FML1PtSmearer.h.

Referenced by FML1PtSmearer(), and smear().

const RandomEngine* FML1PtSmearer::random
private

Definition at line 31 of file FML1PtSmearer.h.

Referenced by smear().

float FML1PtSmearer::resolution[DIMRES]
private

Definition at line 38 of file FML1PtSmearer.h.

Referenced by FML1PtSmearer(), and smear().