Energy scale shifting and smearing module. More...
#include <PhysicsTools/PatAlgos/interface/ObjectEnergyScale.h>
Public Member Functions | |
ObjectEnergyScale (const edm::ParameterSet &iConfig) | |
~ObjectEnergyScale () | |
Private Member Functions | |
float | getSmearing (T &object) |
virtual void | produce (edm::Event &iEvent, const edm::EventSetup &iSetup) |
void | setScale (T &object) |
Private Attributes | |
float | factor_ |
CLHEP::RandGaussQ * | gaussian_ |
float | iniRes_ |
edm::InputTag | objects_ |
float | shiftFactor_ |
bool | useDefaultIniRes_ |
bool | useFixedMass_ |
bool | useIniResByFraction_ |
bool | useWorsenResByFactor_ |
float | worsenRes_ |
Energy scale shifting and smearing module.
This class provides energy scale shifting & smearing to objects with resolutions for systematic error studies. A detailed documentation is found in PhysicsTools/PatAlgos/data/ObjectEnergyScale.cfi
Definition at line 37 of file ObjectEnergyScale.h.
pat::ObjectEnergyScale< T >::ObjectEnergyScale | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 70 of file ObjectEnergyScale.h.
References edm::ParameterSet::getParameter().
{ objects_ = iConfig.getParameter<edm::InputTag>("scaledObject"); useFixedMass_ = iConfig.getParameter<bool> ("fixMass"); shiftFactor_ = iConfig.getParameter<double> ("shiftFactor"); useDefaultIniRes_ = iConfig.getParameter<bool> ("useDefaultInitialResolution"); iniRes_ = iConfig.getParameter<double> ("initialResolution"); useIniResByFraction_ = iConfig.getParameter<bool> ("initialResolutionByFraction"); worsenRes_ = iConfig.getParameter<double> ("worsenResolution"); useWorsenResByFactor_ = iConfig.getParameter<bool> ("worsenResolutionByFactor"); edm::Service<edm::RandomNumberGenerator> rng; CLHEP::HepRandomEngine& engine = rng->getEngine(); gaussian_ = new CLHEP::RandGaussQ(engine); produces<std::vector<T> >(); }
pat::ObjectEnergyScale< T >::~ObjectEnergyScale | ( | ) |
Definition at line 90 of file ObjectEnergyScale.h.
{ delete gaussian_; }
float pat::ObjectEnergyScale< T >::getSmearing | ( | T & | object | ) | [private] |
Returns a smearing factor which is multiplied to the initial value then to get it smeared, sets initial resolution to resolution provided by input object if required and converts the 'worsenResolution' parameter to protect from meaningless final resolution values.
Definition at line 120 of file ObjectEnergyScale.h.
References funct::cos(), relval_parameters_module::energy, max(), funct::pow(), funct::sin(), mathSSE::sqrt(), and theta().
{ // overwrite config file parameter 'initialResolution' if required if ( useDefaultIniRes_ ) { // get initial resolution from input object (and calculate relative initial resolution from absolute value) iniRes_ = (1. / sin(object.theta()) * object.resolutionEt() - object.et() * cos(object.theta()) / std::pow(sin(object.theta()),2) * object.resolutionTheta()) / object.energy(); // conversion of resEt and resTheta into energy resolution } else if ( ! useIniResByFraction_ ) { // calculate relative initial resolution from absolute value iniRes_ = iniRes_ / object.energy(); } // Is 'worsenResolution' a factor or a summand? float finalRes = useWorsenResByFactor_ ? (1.+fabs(1.-fabs(worsenRes_))) * fabs(iniRes_) : fabs(worsenRes_)/object.energy() + fabs(iniRes_); // conversion as protection from "finalRes_<iniRes_" // return smearing factor return std::max( gaussian_->fire(1., sqrt(std::pow(finalRes,2)-std::pow(iniRes_,2))), 0. ); // protection from negative smearing factors }
void pat::ObjectEnergyScale< T >::produce | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [private, virtual] |
Implements edm::EDProducer.
Definition at line 97 of file ObjectEnergyScale.h.
References edm::Event::getByLabel(), i, and edm::Event::put().
{ edm::Handle<std::vector<T> > objectsHandle; iEvent.getByLabel(objects_, objectsHandle); std::vector<T> objects = *objectsHandle; std::auto_ptr<std::vector<T> > objectsVector(new std::vector<T>); objectsVector->reserve(objectsHandle->size()); for ( unsigned int i = 0; i < objects.size(); i++ ) { factor_ = shiftFactor_ * ( objects[i].energy() > 0. ? getSmearing(objects[i]) : 0.); setScale(objects[i]); objectsVector->push_back(objects[i]); } iEvent.put(objectsVector); }
void pat::ObjectEnergyScale< T >::setScale | ( | T & | object | ) | [private] |
Mutliplies the final factor (consisting of shifting and smearing factors) to the object's 4-vector and takes care of preserved masses.
Definition at line 142 of file ObjectEnergyScale.h.
References relval_parameters_module::energy, funct::pow(), and mathSSE::sqrt().
{ if ( factor_ < 0. ) { factor_ = 0.; } // calculate the momentum factor for fixed or not fixed mass float factorMomentum = useFixedMass_ && object.p() > 0. ? sqrt(std::pow(factor_*object.energy(),2)-object.massSqr()) / object.p() : factor_; // set shifted & smeared new 4-vector object.setP4(reco::Particle::LorentzVector(factorMomentum*object.px(), factorMomentum*object.py(), factorMomentum*object.pz(), factor_ *object.energy())); }
float pat::ObjectEnergyScale< T >::factor_ [private] |
Definition at line 52 of file ObjectEnergyScale.h.
CLHEP::RandGaussQ* pat::ObjectEnergyScale< T >::gaussian_ [private] |
Definition at line 61 of file ObjectEnergyScale.h.
float pat::ObjectEnergyScale< T >::iniRes_ [private] |
Definition at line 52 of file ObjectEnergyScale.h.
edm::InputTag pat::ObjectEnergyScale< T >::objects_ [private] |
Definition at line 51 of file ObjectEnergyScale.h.
float pat::ObjectEnergyScale< T >::shiftFactor_ [private] |
Definition at line 52 of file ObjectEnergyScale.h.
bool pat::ObjectEnergyScale< T >::useDefaultIniRes_ [private] |
Definition at line 56 of file ObjectEnergyScale.h.
bool pat::ObjectEnergyScale< T >::useFixedMass_ [private] |
Definition at line 56 of file ObjectEnergyScale.h.
bool pat::ObjectEnergyScale< T >::useIniResByFraction_ [private] |
Definition at line 56 of file ObjectEnergyScale.h.
bool pat::ObjectEnergyScale< T >::useWorsenResByFactor_ [private] |
Definition at line 56 of file ObjectEnergyScale.h.
float pat::ObjectEnergyScale< T >::worsenRes_ [private] |
Definition at line 52 of file ObjectEnergyScale.h.