1 |
dgele |
1.1 |
//
|
2 |
|
|
// $Id: ObjectSpatialResolution.h,v 1.1 2008/03/06 09:23:10 llista Exp $
|
3 |
|
|
//
|
4 |
|
|
|
5 |
|
|
#ifndef PhysicsTools_PatAlgos_ObjectSpatialResolution_h
|
6 |
|
|
#define PhysicsTools_PatAlgos_ObjectSpatialResolution_h
|
7 |
|
|
|
8 |
|
|
/**
|
9 |
|
|
\class pat::ObjectSpatialResolution ObjectSpatialResolution.h "PhysicsTools/PatAlgos/interface/ObjectSpatialResolution.h"
|
10 |
|
|
\brief Energy scale shifting and smearing module
|
11 |
|
|
|
12 |
|
|
This class provides angular smearing to objects with resolutions for
|
13 |
|
|
systematic error studies.
|
14 |
|
|
A detailed documentation is found in
|
15 |
|
|
PhysicsTools/PatAlgos/data/ObjectSpatialResolution.cfi
|
16 |
|
|
|
17 |
|
|
\author Volker Adler
|
18 |
|
|
\version $Id: ObjectSpatialResolution.h,v 1.1 2008/03/06 09:23:10 llista Exp $
|
19 |
|
|
*/
|
20 |
|
|
|
21 |
|
|
|
22 |
|
|
#include <cmath>
|
23 |
|
|
|
24 |
|
|
#include "FWCore/Framework/interface/Frameworkfwd.h"
|
25 |
|
|
#include "FWCore/Framework/interface/EDProducer.h"
|
26 |
|
|
#include "FWCore/Framework/interface/Event.h"
|
27 |
|
|
#include "FWCore/ParameterSet/interface/ParameterSet.h"
|
28 |
|
|
#include "FWCore/ParameterSet/interface/InputTag.h"
|
29 |
|
|
|
30 |
|
|
#include "FWCore/ServiceRegistry/interface/Service.h"
|
31 |
|
|
#include "FWCore/Utilities/interface/RandomNumberGenerator.h"
|
32 |
|
|
#include "CLHEP/Random/RandGaussQ.h"
|
33 |
|
|
#include "DataFormats/Math/interface/normalizedPhi.h"
|
34 |
|
|
|
35 |
|
|
#include <cmath>
|
36 |
|
|
|
37 |
|
|
|
38 |
|
|
namespace pat {
|
39 |
|
|
|
40 |
|
|
|
41 |
|
|
template<class T>
|
42 |
|
|
class ObjectSpatialResolution : public edm::EDProducer {
|
43 |
|
|
|
44 |
|
|
public:
|
45 |
|
|
|
46 |
|
|
explicit ObjectSpatialResolution(const edm::ParameterSet& iConfig);
|
47 |
|
|
~ObjectSpatialResolution();
|
48 |
|
|
|
49 |
|
|
private:
|
50 |
|
|
|
51 |
|
|
virtual void produce(edm::Event& iEvent, const edm::EventSetup& iSetup);
|
52 |
|
|
|
53 |
|
|
void smearAngles(T& object);
|
54 |
|
|
|
55 |
|
|
edm::InputTag objects_;
|
56 |
|
|
double iniResPolar_,
|
57 |
|
|
worsenResPolar_,
|
58 |
|
|
iniResPhi_,
|
59 |
|
|
worsenResPhi_;
|
60 |
|
|
bool useDefaultIniRes_,
|
61 |
|
|
usePolarTheta_,
|
62 |
|
|
useWorsenResPolarByFactor_,
|
63 |
|
|
useWorsenResPhiByFactor_;
|
64 |
|
|
|
65 |
|
|
CLHEP::RandGaussQ* gaussian_;
|
66 |
|
|
|
67 |
|
|
};
|
68 |
|
|
|
69 |
|
|
|
70 |
|
|
}
|
71 |
|
|
|
72 |
|
|
|
73 |
|
|
template<class T>
|
74 |
|
|
pat::ObjectSpatialResolution<T>::ObjectSpatialResolution(const edm::ParameterSet& iConfig)
|
75 |
|
|
{
|
76 |
|
|
objects_ = iConfig.getParameter<edm::InputTag>("movedObject");
|
77 |
|
|
useDefaultIniRes_ = iConfig.getParameter<bool> ("useDefaultInitialResolutions");
|
78 |
|
|
iniResPhi_ = iConfig.getParameter<double> ("initialResolutionPhi");
|
79 |
|
|
worsenResPhi_ = iConfig.getParameter<double> ("worsenResolutionPhi");
|
80 |
|
|
useWorsenResPhiByFactor_ = iConfig.getParameter<bool> ("worsenResolutionPhiByFactor");
|
81 |
|
|
usePolarTheta_ = iConfig.getParameter<bool> ("usePolarTheta");
|
82 |
|
|
iniResPolar_ = iConfig.getParameter<double> ("initialResolutionPolar");
|
83 |
|
|
worsenResPolar_ = iConfig.getParameter<double> ("worsenResolutionPolar");
|
84 |
|
|
useWorsenResPolarByFactor_ = iConfig.getParameter<bool> ("worsenResolutionPolarByFactor");
|
85 |
|
|
|
86 |
|
|
edm::Service<edm::RandomNumberGenerator> rng;
|
87 |
|
|
CLHEP::HepRandomEngine& engine = rng->getEngine();
|
88 |
|
|
gaussian_ = new CLHEP::RandGaussQ(engine);
|
89 |
|
|
|
90 |
|
|
produces<std::vector<T> >();
|
91 |
|
|
}
|
92 |
|
|
|
93 |
|
|
|
94 |
|
|
template<class T>
|
95 |
|
|
pat::ObjectSpatialResolution<T>::~ObjectSpatialResolution()
|
96 |
|
|
{
|
97 |
|
|
delete gaussian_;
|
98 |
|
|
}
|
99 |
|
|
|
100 |
|
|
|
101 |
|
|
template<class T>
|
102 |
|
|
void pat::ObjectSpatialResolution<T>::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
|
103 |
|
|
{
|
104 |
|
|
edm::Handle<std::vector<T> > objectsHandle;
|
105 |
|
|
iEvent.getByLabel(objects_, objectsHandle);
|
106 |
|
|
std::vector<T> objects = *objectsHandle;
|
107 |
|
|
std::auto_ptr<std::vector<T> > objectsVector(new std::vector<T>);
|
108 |
|
|
objectsVector->reserve(objectsHandle->size());
|
109 |
|
|
|
110 |
|
|
for ( unsigned int i = 0; i < objects.size(); i++ ) {
|
111 |
|
|
smearAngles(objects[i]);
|
112 |
|
|
objectsVector->push_back(objects[i]);
|
113 |
|
|
}
|
114 |
|
|
iEvent.put(objectsVector);
|
115 |
|
|
}
|
116 |
|
|
|
117 |
|
|
|
118 |
|
|
/// Sets initial resolution to resolution provided by input object if required,
|
119 |
|
|
/// smears eta/theta and phi and sets the 4-vector accordingl.
|
120 |
|
|
template<class T>
|
121 |
|
|
void pat::ObjectSpatialResolution<T>::smearAngles(T& object)
|
122 |
|
|
{
|
123 |
|
|
// overwrite config file parameters 'initialResolution...' if required
|
124 |
|
|
if ( useDefaultIniRes_ ) {
|
125 |
|
|
// get initial resolutions from input object
|
126 |
|
|
iniResPhi_ = object.resolutionPhi(); // overwrites config file parameter "initialResolution"
|
127 |
|
|
iniResPolar_ = usePolarTheta_ ?
|
128 |
|
|
object.resolutionTheta():
|
129 |
|
|
object.resolutionEta(); // overwrites config file parameter "initialResolution"
|
130 |
|
|
}
|
131 |
|
|
// smear phi
|
132 |
|
|
double finalResPhi = useWorsenResPhiByFactor_ ?
|
133 |
|
|
(1.+fabs(1.-fabs(worsenResPhi_))) * fabs(iniResPhi_):
|
134 |
|
|
fabs(worsenResPhi_) + fabs(iniResPhi_); // conversions as protection from "finalRes_<iniRes_"
|
135 |
|
|
double smearedPhi = normalizedPhi( gaussian_->fire(object.phi(), sqrt(pow(finalResPhi,2)-pow(iniResPhi_,2))) );
|
136 |
|
|
double finalResPolar = useWorsenResPolarByFactor_ ?
|
137 |
|
|
(1.+fabs(1.-fabs(worsenResPolar_))) * fabs(iniResPolar_):
|
138 |
|
|
fabs(worsenResPolar_) + fabs(iniResPolar_); // conversions as protection from "finalRes_<iniRes_"
|
139 |
|
|
// smear theta/eta
|
140 |
|
|
const double thetaMin = 2.*atan(exp(-ROOT::Math::etaMax<double>())); // to be on the safe side; however, etaMax=22765 ==> thetaMin=0 within double precision
|
141 |
|
|
double smearedTheta,
|
142 |
|
|
smearedEta;
|
143 |
|
|
if ( usePolarTheta_ ) {
|
144 |
|
|
smearedTheta = gaussian_->fire(object.theta(), sqrt(pow(finalResPolar,2)-pow(iniResPolar_,2)));
|
145 |
|
|
// 0<theta<Pi needs to be assured to have proper calculation of eta
|
146 |
|
|
while ( fabs(smearedTheta) > M_PI ) smearedTheta = smearedTheta < 0. ?
|
147 |
|
|
smearedTheta + 2.*M_PI:
|
148 |
|
|
smearedTheta - 2.*M_PI;
|
149 |
|
|
if ( smearedTheta < 0. ) {
|
150 |
|
|
smearedTheta = -smearedTheta;
|
151 |
|
|
smearedPhi = normalizedPhi(smearedPhi+M_PI);
|
152 |
|
|
}
|
153 |
|
|
smearedEta = smearedTheta < thetaMin ?
|
154 |
|
|
ROOT::Math::etaMax<double>() :
|
155 |
|
|
( smearedTheta > M_PI-thetaMin ?
|
156 |
|
|
-ROOT::Math::etaMax<double>():
|
157 |
|
|
-log(tan(smearedTheta/2.)) ); // eta not calculable for theta=0,Pi, which could occur
|
158 |
|
|
} else {
|
159 |
|
|
smearedEta = gaussian_->fire(object.eta(), sqrt(pow(finalResPolar,2)-pow(iniResPolar_,2)));
|
160 |
|
|
if ( fabs(smearedEta) > ROOT::Math::etaMax<double>() ) smearedEta = smearedEta < 0. ?
|
161 |
|
|
-ROOT::Math::etaMax<double>():
|
162 |
|
|
ROOT::Math::etaMax<double>();
|
163 |
|
|
smearedTheta = 2. * atan(exp(-smearedEta)); // since exp(x)>0 && atan() returns solution closest to 0, 0<theta<Pi should be assured.
|
164 |
|
|
}
|
165 |
|
|
// set smeared new 4-vector
|
166 |
|
|
math::PtEtaPhiELorentzVector newLorentzVector(object.p()*sin(smearedTheta),
|
167 |
|
|
smearedEta,
|
168 |
|
|
smearedPhi,
|
169 |
|
|
object.energy());
|
170 |
|
|
object.setP4(reco::Particle::LorentzVector(newLorentzVector.Px(),
|
171 |
|
|
newLorentzVector.Py(),
|
172 |
|
|
newLorentzVector.Pz(),
|
173 |
|
|
newLorentzVector.E() ));
|
174 |
|
|
}
|
175 |
|
|
|
176 |
|
|
|
177 |
|
|
#endif
|