Main Page   Namespace List   Class Hierarchy   Compound List   File List   Namespace Members   Compound Members   File Members  

Atlfast::MuonSmearer Class Reference

Provides parameterised muon smearing. More...

#include <MuonSmearer.h>

Inheritance diagram for Atlfast::MuonSmearer:

Inheritance graph
[legend]
Collaboration diagram for Atlfast::MuonSmearer:

Collaboration graph
[legend]
List of all members.

Public Methods

 MuonSmearer (const int aseed, const int lumi, const int keymuo)
virtual ~MuonSmearer ()
virtual HepLorentzVector smear (const HepLorentzVector &avec)
 Smear method decleration provided by ISmearer interface.


Private Methods

std::pair< float, float > resolms (float pt, float eta, float phi)
float resolid1 (float pt, float eta)
std::pair< float, float > resolid2 (float pt, float eta)
float resol1 (float sigmams, float sigmamuon, float sigmaid, float sigmatrack)
float resol2 (float sigmamuon, float sigmatrack)
float resol3 ()
float resol4 ()

Private Attributes

int m_lumi
int m_keymuo
float m_magic1
float m_magic2
float m_magic3
float m_magic4
float m_percent

Detailed Description

Provides parameterised muon smearing.

MuonSmearer honours the ISmearer interface and inherets implementation from DefaultSmearer.

Definition at line 35 of file MuonSmearer.h.


Constructor & Destructor Documentation

Atlfast::MuonSmearer::MuonSmearer const int    aseed,
const int    lumi,
const int    keymuo
[inline]
 

Definition at line 42 of file MuonSmearer.h.

References m_keymuo, m_lumi, m_magic1, m_magic2, m_magic3, m_magic4, and m_percent.

00042                                                                     : 
00043         ISmearer(),       DefaultSmearer(aseed), m_lumi(lumi),    m_keymuo(keymuo), 
00044         m_magic1(0.0005), m_magic2(10./7000.),   m_magic3(0.012), m_magic4(0.02),
00045         m_percent(0.01) { }

virtual Atlfast::MuonSmearer::~MuonSmearer   [inline, virtual]
 

Definition at line 46 of file MuonSmearer.h.

00046 { };

Member Function Documentation

HepLorentzVector Atlfast::MuonSmearer::smear const HepLorentzVector &    avec [virtual]
 

Smear method decleration provided by ISmearer interface.

Reimplemented from Atlfast::DefaultSmearer.

Definition at line 31 of file MuonSmearer.cxx.

References resol1(), resol2(), resol3(), resol4(), resolid1(), resolid2(), and resolms().

00031                                                                    {
00032     // ATLFAST fortran routines   
00033     //        RESMUO(KEYMUO,KEYFUN,PT,ETA,PHI)
00034     //        RESOMU (this had been reduced to merely converting radians to degrees)
00035     // have been converted to C++ here, after restructuring the code.
00036     //        RESOLUMU (in Atlfast/mu_atlfast/resmuo.F) is left in fortran.
00037 
00038     float eta      = avec.pseudoRapidity();
00039     float pt       = avec.perp();
00040     float phi      = avec.phi()*deg;
00041     float sigmams;
00042     float sigmamuon;
00043     float sigmaid;
00044     float sigmatrack;
00045     float sigma=0.;
00046     pair<float, float> sigmapr;
00047     switch (m_keymuo) {
00048     case 1:
00049       while( ((sigmapr=resolms(pt, eta, phi)).first) <-1.) {}
00050       sigma=sigmapr.first;
00051       break;
00052     case 2:
00053       while( (sigma=resolid1(pt, eta))               <-1.) {}
00054       break;
00055     case 3:
00056       while( ((sigmapr=resolms(pt, eta, phi)).first) <-1.) {}
00057       sigmams    = sigmapr.first;
00058       sigmamuon  = sigmapr.second;
00059       while( ((sigmapr=resolid2(pt, eta)).first)     <-1.) {}
00060       sigmaid    = sigmapr.first;
00061       sigmatrack = sigmapr.second;
00062 
00063       sigma=resol1(sigmams, sigmamuon, sigmaid, sigmatrack);
00064       break;
00065     case 4:
00066       while( ((sigmapr=resolms(pt, eta, phi)).first) <-1.) {}
00067       sigmams    = sigmapr.first;
00068       sigmamuon  = sigmapr.second;
00069       while( ((sigmapr=resolid2(pt, eta)).first)     <-1.) {}
00070       sigmaid    = sigmapr.first;
00071       sigmatrack = sigmapr.second;
00072 
00073       sigma=resol2(sigmamuon, sigmatrack);
00074       break;
00075     case 5:
00076       sigma=resol3();
00077       break;
00078     case 6:
00079       sigma=resol4();
00080       
00081     }
00082 
00083     HepLorentzVector bvec(avec.px()/(1.0+sigma), avec.py()/(1.0+sigma), 
00084                           avec.pz()/(1.0+sigma), avec.e() /(1.0+sigma));
00085     return bvec;
00086   }

pair< float, float > Atlfast::MuonSmearer::resolms float    pt,
float    eta,
float    phi
[private]
 

Definition at line 88 of file MuonSmearer.cxx.

References m_percent, Atlfast::DefaultSmearer::randGauss(), and resolumu_().

Referenced by smear().

00088                                                                          {
00089     float cpt=pt;
00090     float ceta=eta;
00091     float cphi=phi;
00092     float resol;
00093     resolumu_(&cpt, &ceta, &cphi, &resol);
00094     double aa=randGauss()->fire();
00095     float sigmamuon    = m_percent*resol;
00096     float sigmams      = aa*sigmamuon;
00097     pair<float, float>   sigmapr(sigmams, sigmamuon);
00098 
00099     return sigmapr;
00100   }

float Atlfast::MuonSmearer::resolid1 float    pt,
float    eta
[private]
 

Definition at line 102 of file MuonSmearer.cxx.

References m_magic1, m_magic2, m_magic3, and Atlfast::DefaultSmearer::randGauss().

Referenced by smear().

00102                                                   {
00103     double aa=randGauss()->fire();
00104     double bb=randGauss()->fire();
00105     float factor = 1. + pow(fabs(eta),m_magic2);
00106     float atrack = m_magic1*factor;
00107     float sigma = atrack*pt*aa + m_magic3*bb;
00108     sigma = sigma*factor;
00109       
00110     return sigma;
00111   }

pair< float, float > Atlfast::MuonSmearer::resolid2 float    pt,
float    eta
[private]
 

Definition at line 113 of file MuonSmearer.cxx.

References m_magic1, m_magic2, m_magic3, and Atlfast::DefaultSmearer::randGauss().

Referenced by smear().

00113                                                                {
00114     double aa=randGauss()->fire();
00115     double bb=randGauss()->fire();
00116     double factor = 1. + pow(fabs(eta),m_magic2);
00117     float atrack = m_magic1*factor;
00118     float sigmatr    = atrack*pt*aa + m_magic3*bb;
00119     float sigmatrack = sqrt( (atrack*pt)*(atrack*pt) + m_magic3*m_magic3);
00120     pair<float, float> sigmapr(sigmatr, sigmatrack);
00121       
00122     return sigmapr;
00123   }

float Atlfast::MuonSmearer::resol1 float    sigmams,
float    sigmamuon,
float    sigmaid,
float    sigmatrack
[private]
 

Definition at line 125 of file MuonSmearer.cxx.

Referenced by smear().

00126                                                                {
00127     double wmuon  = 1./(sigmamuon*sigmamuon);
00128     double wtrack = 1./(sigmatrack*sigmatrack);
00129     double wtot   = wmuon+wtrack;
00130     double corr   = (wmuon*(1.0+sigmams)+wtrack*(1.0+sigmaid))/wtot;
00131     float sigma   = corr-1.0;
00132       
00133     return sigma;
00134   }

float Atlfast::MuonSmearer::resol2 float    sigmamuon,
float    sigmatrack
[private]
 

Definition at line 136 of file MuonSmearer.cxx.

References Atlfast::DefaultSmearer::randGauss().

Referenced by smear().

00136                                                               {
00137     double wmuon  = 1./(sigmamuon*sigmamuon);
00138     double wtrack = 1./(sigmatrack*sigmatrack);
00139     double wtot   = wmuon+wtrack;
00140     double aa     = randGauss()->fire();
00141     float sigma   = aa/sqrt(wtot);
00142       
00143     return sigma;
00144   }

float Atlfast::MuonSmearer::resol3   [private]
 

Definition at line 146 of file MuonSmearer.cxx.

References m_magic4, and Atlfast::DefaultSmearer::randGauss().

Referenced by smear().

00146                              {
00147     float sigma  = m_magic4*randGauss()->fire();
00148       
00149     return sigma;
00150   }

float Atlfast::MuonSmearer::resol4   [private]
 

Definition at line 152 of file MuonSmearer.cxx.

References m_magic4, and Atlfast::DefaultSmearer::randGauss().

Referenced by smear().

00152                              {
00153     float sigmams = m_magic4;
00154     float sigmaid = m_magic4;
00155     float wmuon   = 1.0/(sigmams*sigmams);
00156     float wtrak   = 1.0/(sigmaid*sigmaid);
00157     float wtot    = wmuon + wtrak;
00158     double aa     = randGauss()->fire();
00159     float sigma   = aa/sqrt(wtot);
00160       
00161     return sigma;
00162   }

Member Data Documentation

int Atlfast::MuonSmearer::m_lumi [private]
 

Definition at line 74 of file MuonSmearer.h.

Referenced by MuonSmearer().

int Atlfast::MuonSmearer::m_keymuo [private]
 

Definition at line 75 of file MuonSmearer.h.

Referenced by MuonSmearer().

float Atlfast::MuonSmearer::m_magic1 [private]
 

Definition at line 76 of file MuonSmearer.h.

Referenced by MuonSmearer(), resolid1(), and resolid2().

float Atlfast::MuonSmearer::m_magic2 [private]
 

Definition at line 77 of file MuonSmearer.h.

Referenced by MuonSmearer(), resolid1(), and resolid2().

float Atlfast::MuonSmearer::m_magic3 [private]
 

Definition at line 78 of file MuonSmearer.h.

Referenced by MuonSmearer(), resolid1(), and resolid2().

float Atlfast::MuonSmearer::m_magic4 [private]
 

Definition at line 79 of file MuonSmearer.h.

Referenced by MuonSmearer(), resol3(), and resol4().

float Atlfast::MuonSmearer::m_percent [private]
 

Definition at line 80 of file MuonSmearer.h.

Referenced by MuonSmearer(), and resolms().


The documentation for this class was generated from the following files:
Generated on Tue Mar 18 11:18:58 2003 for AtlfastAlgs by doxygen1.3-rc1