CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Public Attributes
MultipleAlgoIterator Class Reference

#include <MultipleAlgoIterator.h>

Inheritance diagram for MultipleAlgoIterator:
PileUpSubtractor

Public Member Functions

virtual void calculatePedestal (std::vector< fastjet::PseudoJet > const &coll)
 
double getEt (const reco::CandidatePtr &in) const
 
double getEta (const reco::CandidatePtr &in) const
 
 MultipleAlgoIterator (const edm::ParameterSet &iConfig)
 
virtual void offsetCorrectJets ()
 
void rescaleRMS (double s)
 
virtual void subtractPedestal (std::vector< fastjet::PseudoJet > &coll)
 
 ~MultipleAlgoIterator ()
 
- Public Member Functions inherited from PileUpSubtractor
virtual void calculateOrphanInput (std::vector< fastjet::PseudoJet > &orphanInput)
 
virtual double getMeanAtTower (const reco::CandidatePtr &in) const
 
virtual double getPileUpAtTower (const reco::CandidatePtr &in) const
 
virtual double getPileUpEnergy (int ijet) const
 
virtual double getSigmaAtTower (const reco::CandidatePtr &in) const
 
int ieta (const reco::CandidatePtr &in) const
 
int iphi (const reco::CandidatePtr &in) const
 
 PileUpSubtractor (const edm::ParameterSet &iConfig)
 
virtual void reset (std::vector< edm::Ptr< reco::Candidate > > &input, std::vector< fastjet::PseudoJet > &towers, std::vector< fastjet::PseudoJet > &output)
 
virtual void setAlgorithm (ClusterSequencePtr &algorithm)
 
virtual void setupGeometryMap (edm::Event &iEvent, const edm::EventSetup &iSetup)
 
 ~PileUpSubtractor ()
 

Public Attributes

bool dropZeroTowers_
 
bool sumRecHits_
 

Additional Inherited Members

- Public Types inherited from PileUpSubtractor
typedef boost::shared_ptr
< fastjet::ActiveAreaSpec > 
ActiveAreaSpecPtr
 
typedef boost::shared_ptr
< fastjet::ClusterSequence > 
ClusterSequencePtr
 
typedef boost::shared_ptr
< fastjet::RangeDefinition > 
RangeDefPtr
 
- Protected Attributes inherited from PileUpSubtractor
std::vector< HcalDetIdallgeomid_
 
bool doAreaFastjet_
 
bool doRhoFastjet_
 
std::map< int, double > emean_
 
std::map< int, double > esigma_
 
ActiveAreaSpecPtr fjActiveArea_
 
ClusterSequencePtr fjClusterSeq_
 
std::vector< fastjet::PseudoJet > * fjInputs_
 
std::vector< fastjet::PseudoJet > * fjJets_
 
std::vector< fastjet::PseudoJet > fjOriginalInputs_
 
RangeDefPtr fjRangeDef_
 
CaloGeometry const * geo_
 
std::map< int, int > geomtowers_
 
int ietamax_
 
int ietamin_
 
std::vector< edm::Ptr
< reco::Candidate > > * 
inputs_
 
std::vector< double > jetOffset_
 
double jetPtMin_
 
double nSigmaPU_
 
std::map< int, int > ntowersWithJets_
 
double radiusPU_
 
bool reRunAlgo_
 

Detailed Description

Definition at line 6 of file MultipleAlgoIterator.h.

Constructor & Destructor Documentation

MultipleAlgoIterator::MultipleAlgoIterator ( const edm::ParameterSet iConfig)
inline

Definition at line 8 of file MultipleAlgoIterator.h.

8  : PileUpSubtractor(iConfig),
9  sumRecHits_(iConfig.getParameter<bool>("sumRecHits")),
10  dropZeroTowers_(iConfig.getUntrackedParameter<bool>("dropZeroTowers",true))
11  {;}
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
PileUpSubtractor(const edm::ParameterSet &iConfig)
MultipleAlgoIterator::~MultipleAlgoIterator ( )
inline

Definition at line 21 of file MultipleAlgoIterator.h.

21 {;}

Member Function Documentation

void MultipleAlgoIterator::calculatePedestal ( std::vector< fastjet::PseudoJet > const &  coll)
virtual

Reimplemented from PileUpSubtractor.

Definition at line 112 of file MultipleAlgoIterator.cc.

References gt, i, LogDebug, nt, edm::second(), and mathSSE::sqrt().

113 {
114  LogDebug("PileUpSubtractor")<<"The subtractor calculating pedestals...\n";
115 
116  map<int,double> emean2;
117  map<int,int> ntowers;
118 
119  int ietaold = -10000;
120  int ieta0 = -100;
121 
122  // Initial values for emean_, emean2, esigma_, ntowers
123 
124  for(int i = ietamin_; i < ietamax_+1; i++)
125  {
126  emean_[i] = 0.;
127  emean2[i] = 0.;
128  esigma_[i] = 0.;
129  ntowers[i] = 0;
130  }
131 
132  for (vector<fastjet::PseudoJet>::const_iterator input_object = coll.begin (),
133  fjInputsEnd = coll.end();
134  input_object != fjInputsEnd; ++input_object) {
135 
136  const reco::CandidatePtr & originalTower=(*inputs_)[ input_object->user_index()];
137  ieta0 = ieta( originalTower );
138  double Original_Et = originalTower->et();
139  if(sumRecHits_){
140  Original_Et = getEt(originalTower);
141  }
142 
143  if( ieta0-ietaold != 0 )
144  {
145  emean_[ieta0] = emean_[ieta0]+Original_Et;
146  emean2[ieta0] = emean2[ieta0]+Original_Et*Original_Et;
147  ntowers[ieta0] = 1;
148  ietaold = ieta0;
149  }
150  else
151  {
152  emean_[ieta0] = emean_[ieta0]+Original_Et;
153  emean2[ieta0] = emean2[ieta0]+Original_Et*Original_Et;
154  ntowers[ieta0]++;
155  }
156 
157  }
158 
159  for(map<int,int>::const_iterator gt = geomtowers_.begin(); gt != geomtowers_.end(); gt++)
160  {
161 
162  int it = (*gt).first;
163 
164  double e1 = (*(emean_.find(it))).second;
165  double e2 = (*emean2.find(it)).second;
166  int nt = (*gt).second - (*(ntowersWithJets_.find(it))).second;
167 
168  LogDebug("PileUpSubtractor")<<" ieta : "<<it<<" number of towers : "<<nt<<" e1 : "<<e1<<" e2 : "<<e2<<"\n";
169 
170  if(nt > 0) {
171  emean_[it] = e1/nt;
172  double eee = e2/nt - e1*e1/(nt*nt);
173  if(eee<0.) eee = 0.;
174  esigma_[it] = nSigmaPU_*sqrt(eee);
175  }
176  else
177  {
178  emean_[it] = 0.;
179  esigma_[it] = 0.;
180  }
181  LogDebug("PileUpSubtractor")<<" ieta : "<<it<<" Pedestals : "<<emean_[it]<<" "<<esigma_[it]<<"\n";
182  }
183 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
double getEt(const reco::CandidatePtr &in) const
std::map< int, double > esigma_
std::map< int, int > geomtowers_
int ieta(const reco::CandidatePtr &in) const
U second(std::pair< T, U > const &p)
std::map< int, int > ntowersWithJets_
T sqrt(T t)
Definition: SSEVec.h:28
int nt
Definition: AMPTWrapper.h:32
JetCorrectorParametersCollection coll
Definition: classes.h:14
std::map< int, double > emean_
DTDigi & gt
double MultipleAlgoIterator::getEt ( const reco::CandidatePtr in) const

Definition at line 189 of file MultipleAlgoIterator.cc.

References CaloTower::emEnergy(), relval_parameters_module::energy, edm::Ptr< T >::get(), CaloTower::hadEnergy(), CaloTower::id(), pos, funct::sin(), and PV3DBase< T, PVType, FrameType >::theta().

189  {
190  const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
191  const GlobalPoint& pos=geo_->getPosition(ctc->id());
192  double energy = ctc->emEnergy() + ctc->hadEnergy();
193  double et = energy*sin(pos.theta());
194  return et;
195 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
Definition: PV3DBase.h:69
T const * get() const
Returns C++ pointer to the item.
Definition: Ptr.h:148
double emEnergy() const
Definition: CaloTower.h:79
const GlobalPoint & getPosition(const DetId &id) const
Get the position of a given detector id.
Definition: CaloGeometry.cc:68
double hadEnergy() const
Definition: CaloTower.h:80
CaloTowerDetId id() const
Definition: CaloTower.h:72
CaloGeometry const * geo_
double MultipleAlgoIterator::getEta ( const reco::CandidatePtr in) const

Definition at line 197 of file MultipleAlgoIterator.cc.

References PV3DBase< T, PVType, FrameType >::eta(), eta(), edm::Ptr< T >::get(), CaloTower::id(), and pos.

197  {
198  const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
199  const GlobalPoint& pos=geo_->getPosition(ctc->id());
200  double eta = pos.eta();
201  return eta;
202 }
T eta() const
T const * get() const
Returns C++ pointer to the item.
Definition: Ptr.h:148
const GlobalPoint & getPosition(const DetId &id) const
Get the position of a given detector id.
Definition: CaloGeometry.cc:68
CaloTowerDetId id() const
Definition: CaloTower.h:72
T eta() const
Definition: PV3DBase.h:70
CaloGeometry const * geo_
void MultipleAlgoIterator::offsetCorrectJets ( )
virtual

Reimplemented from PileUpSubtractor.

Definition at line 16 of file MultipleAlgoIterator.cc.

References def, LogDebug, dt_dqm_sourceclient_common_cff::reco, and edm::second().

17 {
18 
19  LogDebug("PileUpSubtractor")<<"The subtractor correcting jets...\n";
20  jetOffset_.clear();
21 
22  using namespace reco;
23 
24  (*fjInputs_) = fjOriginalInputs_;
27  const fastjet::JetDefinition& def = fjClusterSeq_->jet_def();
28  if ( !doAreaFastjet_ && !doRhoFastjet_) {
29  fastjet::ClusterSequence newseq( *fjInputs_, def );
30  (*fjClusterSeq_) = newseq;
31  } else {
32  fastjet::ClusterSequenceArea newseq( *fjInputs_, def , *fjActiveArea_ );
33  (*fjClusterSeq_) = newseq;
34  }
35 
36  (*fjJets_) = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
37 
38  jetOffset_.reserve(fjJets_->size());
39 
40  vector<fastjet::PseudoJet>::iterator pseudojetTMP = fjJets_->begin (),
41  jetsEnd = fjJets_->end();
42  for (; pseudojetTMP != jetsEnd; ++pseudojetTMP) {
43 
44  int ijet = pseudojetTMP - fjJets_->begin();
45  jetOffset_[ijet] = 0;
46 
47  std::vector<fastjet::PseudoJet> towers =
48  sorted_by_pt(fjClusterSeq_->constituents(*pseudojetTMP));
49 
50  double newjetet = 0.;
51  for(vector<fastjet::PseudoJet>::const_iterator ito = towers.begin(),
52  towEnd = towers.end();
53  ito != towEnd;
54  ++ito)
55  {
56  const reco::CandidatePtr& originalTower = (*inputs_)[ito->user_index()];
57  int it = ieta( originalTower );
58  double Original_Et = originalTower->et();
59  double etnew = Original_Et - (*emean_.find(it)).second - (*esigma_.find(it)).second;
60  if(etnew < 0.) etnew = 0;
61  newjetet = newjetet + etnew;
62  jetOffset_[ijet] += Original_Et - etnew;
63  }
64  }
65 }
#define LogDebug(id)
std::vector< double > jetOffset_
std::vector< fastjet::PseudoJet > * fjJets_
std::map< int, double > esigma_
std::vector< fastjet::PseudoJet > fjOriginalInputs_
virtual void subtractPedestal(std::vector< fastjet::PseudoJet > &coll)
int ieta(const reco::CandidatePtr &in) const
U second(std::pair< T, U > const &p)
std::vector< fastjet::PseudoJet > * fjInputs_
ActiveAreaSpecPtr fjActiveArea_
ClusterSequencePtr fjClusterSeq_
std::map< int, double > emean_
JetCorrectorParameters::Definitions def
Definition: classes.h:10
void MultipleAlgoIterator::rescaleRMS ( double  s)

Definition at line 8 of file MultipleAlgoIterator.cc.

8  {
9  for ( std::map<int, double>::iterator iter = esigma_.begin();
10  iter != esigma_.end(); ++iter ){
11  iter->second = s*(iter->second);
12  }
13 }
std::map< int, double > esigma_
string s
Definition: asciidump.py:422
void MultipleAlgoIterator::subtractPedestal ( std::vector< fastjet::PseudoJet > &  coll)
virtual

Reimplemented from PileUpSubtractor.

Definition at line 67 of file MultipleAlgoIterator.cc.

References getHLTprescales::index, and LogDebug.

68 {
69 
70  LogDebug("PileUpSubtractor")<<"The subtractor subtracting pedestals...\n";
71 
72  int it = -100;
73  int ip = -100;
74 
75  vector<fastjet::PseudoJet> newcoll;
76 
77  for (vector<fastjet::PseudoJet>::iterator input_object = coll.begin (),
78  fjInputsEnd = coll.end();
79  input_object != fjInputsEnd; ++input_object) {
80 
81  reco::CandidatePtr const & itow = (*inputs_)[ input_object->user_index() ];
82 
83  it = ieta( itow );
84  ip = iphi( itow );
85 
86  double Original_Et = itow->et();
87  if(sumRecHits_){
88  Original_Et = getEt(itow);
89  }
90 
91  double etnew = Original_Et - (*(emean_.find(it))).second - (*(esigma_.find(it))).second;
92  float mScale = etnew/input_object->Et();
93  if(etnew < 0.) mScale = 0.;
94 
95  math::XYZTLorentzVectorD towP4(input_object->px()*mScale, input_object->py()*mScale,
96  input_object->pz()*mScale, input_object->e()*mScale);
97 
98  int index = input_object->user_index();
99  input_object->reset ( towP4.px(),
100  towP4.py(),
101  towP4.pz(),
102  towP4.energy() );
103  input_object->set_user_index(index);
104 
105  if(etnew > 0. && dropZeroTowers_) newcoll.push_back(*input_object);
106  }
107 
108  if(dropZeroTowers_) coll = newcoll;
109 
110 }
#define LogDebug(id)
double getEt(const reco::CandidatePtr &in) const
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:15
std::map< int, double > esigma_
int ieta(const reco::CandidatePtr &in) const
int iphi(const reco::CandidatePtr &in) const
JetCorrectorParametersCollection coll
Definition: classes.h:14
std::map< int, double > emean_

Member Data Documentation

bool MultipleAlgoIterator::dropZeroTowers_

Definition at line 20 of file MultipleAlgoIterator.h.

bool MultipleAlgoIterator::sumRecHits_

Definition at line 19 of file MultipleAlgoIterator.h.