Go to the documentation of this file.00001
00002
00003
00004
00005 #ifndef PhysicsTools_PatAlgos_ObjectEnergyScale_h
00006 #define PhysicsTools_PatAlgos_ObjectEnergyScale_h
00007
00022 #include "FWCore/Framework/interface/Frameworkfwd.h"
00023 #include "FWCore/Framework/interface/EDProducer.h"
00024 #include "FWCore/Framework/interface/Event.h"
00025 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00026 #include "FWCore/Utilities/interface/InputTag.h"
00027
00028 #include "FWCore/ServiceRegistry/interface/Service.h"
00029 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
00030 #include "CLHEP/Random/RandGaussQ.h"
00031
00032
00033 namespace pat {
00034
00035
00036 template<class T>
00037 class ObjectEnergyScale : public edm::EDProducer {
00038
00039 public:
00040
00041 explicit ObjectEnergyScale(const edm::ParameterSet& iConfig);
00042 ~ObjectEnergyScale();
00043
00044 private:
00045
00046 virtual void produce( edm::Event& iEvent, const edm::EventSetup& iSetup);
00047
00048 float getSmearing(T& object);
00049 void setScale(T& object);
00050
00051 edm::InputTag objects_;
00052 float factor_,
00053 shiftFactor_,
00054 iniRes_,
00055 worsenRes_;
00056 bool useFixedMass_,
00057 useDefaultIniRes_,
00058 useIniResByFraction_,
00059 useWorsenResByFactor_;
00060
00061 CLHEP::RandGaussQ* gaussian_;
00062
00063 };
00064
00065
00066 }
00067
00068
00069 template<class T>
00070 pat::ObjectEnergyScale<T>::ObjectEnergyScale(const edm::ParameterSet& iConfig)
00071 {
00072 objects_ = iConfig.getParameter<edm::InputTag>("scaledObject");
00073 useFixedMass_ = iConfig.getParameter<bool> ("fixMass");
00074 shiftFactor_ = iConfig.getParameter<double> ("shiftFactor");
00075 useDefaultIniRes_ = iConfig.getParameter<bool> ("useDefaultInitialResolution");
00076 iniRes_ = iConfig.getParameter<double> ("initialResolution");
00077 useIniResByFraction_ = iConfig.getParameter<bool> ("initialResolutionByFraction");
00078 worsenRes_ = iConfig.getParameter<double> ("worsenResolution");
00079 useWorsenResByFactor_ = iConfig.getParameter<bool> ("worsenResolutionByFactor");
00080
00081 edm::Service<edm::RandomNumberGenerator> rng;
00082 CLHEP::HepRandomEngine& engine = rng->getEngine();
00083 gaussian_ = new CLHEP::RandGaussQ(engine);
00084
00085 produces<std::vector<T> >();
00086 }
00087
00088
00089 template<class T>
00090 pat::ObjectEnergyScale<T>::~ObjectEnergyScale()
00091 {
00092 delete gaussian_;
00093 }
00094
00095
00096 template<class T>
00097 void pat::ObjectEnergyScale<T>::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00098 {
00099 edm::Handle<std::vector<T> > objectsHandle;
00100 iEvent.getByLabel(objects_, objectsHandle);
00101 std::vector<T> objects = *objectsHandle;
00102 std::auto_ptr<std::vector<T> > objectsVector(new std::vector<T>);
00103 objectsVector->reserve(objectsHandle->size());
00104
00105 for ( unsigned int i = 0; i < objects.size(); i++ ) {
00106 factor_ = shiftFactor_ * ( objects[i].energy() > 0. ?
00107 getSmearing(objects[i]) :
00108 0.);
00109 setScale(objects[i]);
00110 objectsVector->push_back(objects[i]);
00111 }
00112 iEvent.put(objectsVector);
00113 }
00114
00115
00119 template<class T>
00120 float pat::ObjectEnergyScale<T>::getSmearing(T& object)
00121 {
00122
00123 if ( useDefaultIniRes_ ) {
00124
00125 iniRes_ = (1. / sin(object.theta()) * object.resolutionEt() - object.et() * cos(object.theta()) / std::pow(sin(object.theta()),2) * object.resolutionTheta()) / object.energy();
00126 } else if ( ! useIniResByFraction_ ) {
00127
00128 iniRes_ = iniRes_ / object.energy();
00129 }
00130
00131 float finalRes = useWorsenResByFactor_ ?
00132 (1.+fabs(1.-fabs(worsenRes_))) * fabs(iniRes_) :
00133 fabs(worsenRes_)/object.energy() + fabs(iniRes_);
00134
00135 return std::max( gaussian_->fire(1., sqrt(std::pow(finalRes,2)-std::pow(iniRes_,2))), 0. );
00136 }
00137
00138
00141 template<class T>
00142 void pat::ObjectEnergyScale<T>::setScale(T& object)
00143 {
00144 if ( factor_ < 0. ) {
00145 factor_ = 0.;
00146 }
00147
00148 float factorMomentum = useFixedMass_ && object.p() > 0. ?
00149 sqrt(std::pow(factor_*object.energy(),2)-object.massSqr()) / object.p() :
00150 factor_;
00151
00152 object.setP4(reco::Particle::LorentzVector(factorMomentum*object.px(),
00153 factorMomentum*object.py(),
00154 factorMomentum*object.pz(),
00155 factor_ *object.energy()));
00156 }
00157
00158
00159 #endif