CMS 3D CMS Logo

pat::ObjectSpatialResolution< T > Class Template Reference

Energy scale shifting and smearing module. More...

#include <PhysicsTools/PatAlgos/interface/ObjectSpatialResolution.h>

Inheritance diagram for pat::ObjectSpatialResolution< T >:

edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

 ObjectSpatialResolution (const edm::ParameterSet &iConfig)
 ~ObjectSpatialResolution ()

Private Member Functions

virtual void produce (edm::Event &iEvent, const edm::EventSetup &iSetup)
void smearAngles (T &object)
 Sets initial resolution to resolution provided by input object if required, smears eta/theta and phi and sets the 4-vector accordingl.

Private Attributes

CLHEP::RandGaussQ * gaussian_
double iniResPhi_
double iniResPolar_
edm::InputTag objects_
bool useDefaultIniRes_
bool usePolarTheta_
bool useWorsenResPhiByFactor_
bool useWorsenResPolarByFactor_
double worsenResPhi_
double worsenResPolar_


Detailed Description

template<class T>
class pat::ObjectSpatialResolution< T >

Energy scale shifting and smearing module.

This class provides angular smearing to objects with resolutions for systematic error studies. A detailed documentation is found in PhysicsTools/PatAlgos/data/ObjectSpatialResolution.cfi

Author:
Volker Adler
Version:
Id
ObjectSpatialResolution.h,v 1.1 2008/03/06 09:23:10 llista Exp

Definition at line 42 of file ObjectSpatialResolution.h.


Constructor & Destructor Documentation

template<class T>
pat::ObjectSpatialResolution< T >::ObjectSpatialResolution ( const edm::ParameterSet iConfig  )  [inline, explicit]

Definition at line 74 of file ObjectSpatialResolution.h.

References pat::ObjectSpatialResolution< T >::gaussian_, edm::ParameterSet::getParameter(), pat::ObjectSpatialResolution< T >::iniResPhi_, pat::ObjectSpatialResolution< T >::iniResPolar_, pat::ObjectSpatialResolution< T >::objects_, pat::ObjectSpatialResolution< T >::useDefaultIniRes_, pat::ObjectSpatialResolution< T >::usePolarTheta_, pat::ObjectSpatialResolution< T >::useWorsenResPhiByFactor_, pat::ObjectSpatialResolution< T >::useWorsenResPolarByFactor_, pat::ObjectSpatialResolution< T >::worsenResPhi_, and pat::ObjectSpatialResolution< T >::worsenResPolar_.

00075 {
00076   objects_                   = iConfig.getParameter<edm::InputTag>("movedObject");
00077   useDefaultIniRes_          = iConfig.getParameter<bool>         ("useDefaultInitialResolutions");
00078   iniResPhi_                 = iConfig.getParameter<double>       ("initialResolutionPhi");
00079   worsenResPhi_              = iConfig.getParameter<double>       ("worsenResolutionPhi");
00080   useWorsenResPhiByFactor_   = iConfig.getParameter<bool>         ("worsenResolutionPhiByFactor");
00081   usePolarTheta_             = iConfig.getParameter<bool>         ("usePolarTheta");
00082   iniResPolar_               = iConfig.getParameter<double>       ("initialResolutionPolar");
00083   worsenResPolar_            = iConfig.getParameter<double>       ("worsenResolutionPolar");
00084   useWorsenResPolarByFactor_ = iConfig.getParameter<bool>         ("worsenResolutionPolarByFactor");
00085 
00086   edm::Service<edm::RandomNumberGenerator> rng;
00087   CLHEP::HepRandomEngine& engine = rng->getEngine();
00088   gaussian_ = new CLHEP::RandGaussQ(engine);
00089 
00090   produces<std::vector<T> >();
00091 }

template<class T>
pat::ObjectSpatialResolution< T >::~ObjectSpatialResolution (  )  [inline]

Definition at line 95 of file ObjectSpatialResolution.h.

References pat::ObjectSpatialResolution< T >::gaussian_.

00096 {
00097   delete gaussian_;
00098 }


Member Function Documentation

template<class T>
void pat::ObjectSpatialResolution< T >::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
) [inline, private, virtual]

Implements edm::EDProducer.

Definition at line 102 of file ObjectSpatialResolution.h.

References edm::Event::getByLabel(), i, pat::ObjectSpatialResolution< T >::objects_, edm::Event::put(), and pat::ObjectSpatialResolution< T >::smearAngles().

00103 {
00104   edm::Handle<std::vector<T> > objectsHandle;
00105   iEvent.getByLabel(objects_, objectsHandle);
00106   std::vector<T> objects = *objectsHandle;
00107   std::auto_ptr<std::vector<T> > objectsVector(new std::vector<T>);
00108   objectsVector->reserve(objectsHandle->size());
00109 
00110   for ( unsigned int i = 0; i < objects.size(); i++ ) {
00111     smearAngles(objects[i]);
00112     objectsVector->push_back(objects[i]);
00113   }
00114   iEvent.put(objectsVector);
00115 }

template<class T>
void pat::ObjectSpatialResolution< T >::smearAngles ( T &  object  )  [inline, private]

Sets initial resolution to resolution provided by input object if required, smears eta/theta and phi and sets the 4-vector accordingl.

Definition at line 121 of file ObjectSpatialResolution.h.

References relval_parameters_module::energy, eta, funct::exp(), pat::ObjectSpatialResolution< T >::gaussian_, pat::ObjectSpatialResolution< T >::iniResPhi_, pat::ObjectSpatialResolution< T >::iniResPolar_, funct::log(), normalizedPhi(), p, funct::pow(), funct::sin(), funct::sqrt(), funct::tan(), theta, pat::ObjectSpatialResolution< T >::useDefaultIniRes_, pat::ObjectSpatialResolution< T >::usePolarTheta_, pat::ObjectSpatialResolution< T >::useWorsenResPhiByFactor_, pat::ObjectSpatialResolution< T >::useWorsenResPolarByFactor_, pat::ObjectSpatialResolution< T >::worsenResPhi_, and pat::ObjectSpatialResolution< T >::worsenResPolar_.

Referenced by pat::ObjectSpatialResolution< T >::produce().

00122 {
00123   // overwrite config file parameters 'initialResolution...' if required
00124   if ( useDefaultIniRes_ ) {
00125     // get initial resolutions from input object
00126     iniResPhi_   = object.resolutionPhi();    // overwrites config file parameter "initialResolution"
00127     iniResPolar_ = usePolarTheta_      ?
00128                    object.resolutionTheta():
00129                    object.resolutionEta();    // overwrites config file parameter "initialResolution"
00130   }
00131   // smear phi
00132   double finalResPhi = useWorsenResPhiByFactor_                            ?
00133                        (1.+fabs(1.-fabs(worsenResPhi_))) * fabs(iniResPhi_):
00134                        fabs(worsenResPhi_) + fabs(iniResPhi_);               // conversions as protection from "finalRes_<iniRes_"
00135   double smearedPhi = normalizedPhi( gaussian_->fire(object.phi(), sqrt(pow(finalResPhi,2)-pow(iniResPhi_,2))) );
00136   double finalResPolar = useWorsenResPolarByFactor_                              ?
00137                          (1.+fabs(1.-fabs(worsenResPolar_))) * fabs(iniResPolar_):
00138                          fabs(worsenResPolar_) + fabs(iniResPolar_);               // conversions as protection from "finalRes_<iniRes_"
00139   // smear theta/eta
00140   const double thetaMin = 2.*atan(exp(-ROOT::Math::etaMax<double>())); // to be on the safe side; however, etaMax=22765 ==> thetaMin=0 within double precision
00141   double smearedTheta,
00142          smearedEta;
00143   if ( usePolarTheta_ ) {
00144     smearedTheta = gaussian_->fire(object.theta(), sqrt(pow(finalResPolar,2)-pow(iniResPolar_,2)));
00145     // 0<theta<Pi needs to be assured to have proper calculation of eta
00146     while ( fabs(smearedTheta) > M_PI ) smearedTheta = smearedTheta < 0.     ?
00147                                                        smearedTheta + 2.*M_PI:
00148                                                        smearedTheta - 2.*M_PI;
00149     if ( smearedTheta < 0. ) {
00150       smearedTheta = -smearedTheta;
00151       smearedPhi   = normalizedPhi(smearedPhi+M_PI);
00152     }
00153     smearedEta = smearedTheta < thetaMin          ?
00154                  ROOT::Math::etaMax<double>()     :
00155                  ( smearedTheta > M_PI-thetaMin ?
00156                    -ROOT::Math::etaMax<double>():
00157                    -log(tan(smearedTheta/2.))    ); // eta not calculable for theta=0,Pi, which could occur
00158   } else {
00159     smearedEta = gaussian_->fire(object.eta(), sqrt(pow(finalResPolar,2)-pow(iniResPolar_,2)));
00160     if ( fabs(smearedEta) > ROOT::Math::etaMax<double>() ) smearedEta = smearedEta < 0.              ?
00161                                                                         -ROOT::Math::etaMax<double>():
00162                                                                          ROOT::Math::etaMax<double>();
00163     smearedTheta = 2. * atan(exp(-smearedEta)); // since exp(x)>0 && atan() returns solution closest to 0, 0<theta<Pi should be assured.
00164   }
00165   // set smeared new 4-vector
00166   math::PtEtaPhiELorentzVector newLorentzVector(object.p()*sin(smearedTheta),
00167                                                 smearedEta,
00168                                                 smearedPhi,
00169                                                 object.energy());
00170   object.setP4(reco::Particle::LorentzVector(newLorentzVector.Px(),
00171                                              newLorentzVector.Py(),
00172                                              newLorentzVector.Pz(),
00173                                              newLorentzVector.E() ));
00174 }


Member Data Documentation

template<class T>
CLHEP::RandGaussQ* pat::ObjectSpatialResolution< T >::gaussian_ [private]

Definition at line 65 of file ObjectSpatialResolution.h.

Referenced by pat::ObjectSpatialResolution< T >::ObjectSpatialResolution(), pat::ObjectSpatialResolution< T >::smearAngles(), and pat::ObjectSpatialResolution< T >::~ObjectSpatialResolution().

template<class T>
double pat::ObjectSpatialResolution< T >::iniResPhi_ [private]

Definition at line 56 of file ObjectSpatialResolution.h.

Referenced by pat::ObjectSpatialResolution< T >::ObjectSpatialResolution(), and pat::ObjectSpatialResolution< T >::smearAngles().

template<class T>
double pat::ObjectSpatialResolution< T >::iniResPolar_ [private]

Definition at line 56 of file ObjectSpatialResolution.h.

Referenced by pat::ObjectSpatialResolution< T >::ObjectSpatialResolution(), and pat::ObjectSpatialResolution< T >::smearAngles().

template<class T>
edm::InputTag pat::ObjectSpatialResolution< T >::objects_ [private]

Definition at line 55 of file ObjectSpatialResolution.h.

Referenced by pat::ObjectSpatialResolution< T >::ObjectSpatialResolution(), and pat::ObjectSpatialResolution< T >::produce().

template<class T>
bool pat::ObjectSpatialResolution< T >::useDefaultIniRes_ [private]

Definition at line 60 of file ObjectSpatialResolution.h.

Referenced by pat::ObjectSpatialResolution< T >::ObjectSpatialResolution(), and pat::ObjectSpatialResolution< T >::smearAngles().

template<class T>
bool pat::ObjectSpatialResolution< T >::usePolarTheta_ [private]

Definition at line 60 of file ObjectSpatialResolution.h.

Referenced by pat::ObjectSpatialResolution< T >::ObjectSpatialResolution(), and pat::ObjectSpatialResolution< T >::smearAngles().

template<class T>
bool pat::ObjectSpatialResolution< T >::useWorsenResPhiByFactor_ [private]

Definition at line 60 of file ObjectSpatialResolution.h.

Referenced by pat::ObjectSpatialResolution< T >::ObjectSpatialResolution(), and pat::ObjectSpatialResolution< T >::smearAngles().

template<class T>
bool pat::ObjectSpatialResolution< T >::useWorsenResPolarByFactor_ [private]

Definition at line 60 of file ObjectSpatialResolution.h.

Referenced by pat::ObjectSpatialResolution< T >::ObjectSpatialResolution(), and pat::ObjectSpatialResolution< T >::smearAngles().

template<class T>
double pat::ObjectSpatialResolution< T >::worsenResPhi_ [private]

Definition at line 56 of file ObjectSpatialResolution.h.

Referenced by pat::ObjectSpatialResolution< T >::ObjectSpatialResolution(), and pat::ObjectSpatialResolution< T >::smearAngles().

template<class T>
double pat::ObjectSpatialResolution< T >::worsenResPolar_ [private]

Definition at line 56 of file ObjectSpatialResolution.h.

Referenced by pat::ObjectSpatialResolution< T >::ObjectSpatialResolution(), and pat::ObjectSpatialResolution< T >::smearAngles().


The documentation for this class was generated from the following file:
Generated on Tue Jun 9 18:49:42 2009 for CMSSW by  doxygen 1.5.4