Go to the documentation of this file.00001 #include "RecoHI/HiJetAlgos/interface/MultipleAlgoIterator.h"
00002 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00003 #include "DataFormats/Candidate/interface/Candidate.h"
00004
00005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00006 using namespace std;
00007
00008 void MultipleAlgoIterator::rescaleRMS(double s){
00009 for ( std::map<int, double>::iterator iter = esigma_.begin();
00010 iter != esigma_.end(); ++iter ){
00011 iter->second = s*(iter->second);
00012 }
00013 }
00014
00015
00016 void MultipleAlgoIterator::offsetCorrectJets()
00017 {
00018
00019 LogDebug("PileUpSubtractor")<<"The subtractor correcting jets...\n";
00020 jetOffset_.clear();
00021
00022 using namespace reco;
00023
00024 (*fjInputs_) = fjOriginalInputs_;
00025 rescaleRMS(nSigmaPU_);
00026 subtractPedestal(*fjInputs_);
00027 const fastjet::JetDefinition& def = fjClusterSeq_->jet_def();
00028 if ( !doAreaFastjet_ && !doRhoFastjet_) {
00029 fastjet::ClusterSequence newseq( *fjInputs_, def );
00030 (*fjClusterSeq_) = newseq;
00031 } else {
00032 fastjet::ClusterSequenceArea newseq( *fjInputs_, def , *fjActiveArea_ );
00033 (*fjClusterSeq_) = newseq;
00034 }
00035
00036 (*fjJets_) = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
00037
00038 jetOffset_.reserve(fjJets_->size());
00039
00040 vector<fastjet::PseudoJet>::iterator pseudojetTMP = fjJets_->begin (),
00041 jetsEnd = fjJets_->end();
00042 for (; pseudojetTMP != jetsEnd; ++pseudojetTMP) {
00043
00044 int ijet = pseudojetTMP - fjJets_->begin();
00045 jetOffset_[ijet] = 0;
00046
00047 std::vector<fastjet::PseudoJet> towers =
00048 sorted_by_pt(fjClusterSeq_->constituents(*pseudojetTMP));
00049
00050 double newjetet = 0.;
00051 for(vector<fastjet::PseudoJet>::const_iterator ito = towers.begin(),
00052 towEnd = towers.end();
00053 ito != towEnd;
00054 ++ito)
00055 {
00056 const reco::CandidatePtr& originalTower = (*inputs_)[ito->user_index()];
00057 int it = ieta( originalTower );
00058 double Original_Et = originalTower->et();
00059 double etnew = Original_Et - (*emean_.find(it)).second - (*esigma_.find(it)).second;
00060 if(etnew < 0.) etnew = 0;
00061 newjetet = newjetet + etnew;
00062 jetOffset_[ijet] += Original_Et - etnew;
00063 }
00064 }
00065 }
00066
00067 void MultipleAlgoIterator::subtractPedestal(vector<fastjet::PseudoJet> & coll)
00068 {
00069
00070 LogDebug("PileUpSubtractor")<<"The subtractor subtracting pedestals...\n";
00071
00072 int it = -100;
00073 int ip = -100;
00074
00075 vector<fastjet::PseudoJet> newcoll;
00076
00077 for (vector<fastjet::PseudoJet>::iterator input_object = coll.begin (),
00078 fjInputsEnd = coll.end();
00079 input_object != fjInputsEnd; ++input_object) {
00080
00081 reco::CandidatePtr const & itow = (*inputs_)[ input_object->user_index() ];
00082
00083 it = ieta( itow );
00084 ip = iphi( itow );
00085
00086 double Original_Et = itow->et();
00087 if(sumRecHits_){
00088 Original_Et = getEt(itow);
00089 }
00090
00091 double etnew = Original_Et - (*(emean_.find(it))).second - (*(esigma_.find(it))).second;
00092 float mScale = etnew/input_object->Et();
00093 if(etnew < 0.) mScale = 0.;
00094
00095 math::XYZTLorentzVectorD towP4(input_object->px()*mScale, input_object->py()*mScale,
00096 input_object->pz()*mScale, input_object->e()*mScale);
00097
00098 int index = input_object->user_index();
00099 input_object->reset ( towP4.px(),
00100 towP4.py(),
00101 towP4.pz(),
00102 towP4.energy() );
00103 input_object->set_user_index(index);
00104
00105 if(etnew > 0. && dropZeroTowers_) newcoll.push_back(*input_object);
00106 }
00107
00108 if(dropZeroTowers_) coll = newcoll;
00109
00110 }
00111
00112 void MultipleAlgoIterator::calculatePedestal( vector<fastjet::PseudoJet> const & coll )
00113 {
00114 LogDebug("PileUpSubtractor")<<"The subtractor calculating pedestals...\n";
00115
00116 map<int,double> emean2;
00117 map<int,int> ntowers;
00118
00119 int ietaold = -10000;
00120 int ieta0 = -100;
00121
00122
00123
00124 for(int i = ietamin_; i < ietamax_+1; i++)
00125 {
00126 emean_[i] = 0.;
00127 emean2[i] = 0.;
00128 esigma_[i] = 0.;
00129 ntowers[i] = 0;
00130 }
00131
00132 for (vector<fastjet::PseudoJet>::const_iterator input_object = coll.begin (),
00133 fjInputsEnd = coll.end();
00134 input_object != fjInputsEnd; ++input_object) {
00135
00136 const reco::CandidatePtr & originalTower=(*inputs_)[ input_object->user_index()];
00137 ieta0 = ieta( originalTower );
00138 double Original_Et = originalTower->et();
00139 if(sumRecHits_){
00140 Original_Et = getEt(originalTower);
00141 }
00142
00143 if( ieta0-ietaold != 0 )
00144 {
00145 emean_[ieta0] = emean_[ieta0]+Original_Et;
00146 emean2[ieta0] = emean2[ieta0]+Original_Et*Original_Et;
00147 ntowers[ieta0] = 1;
00148 ietaold = ieta0;
00149 }
00150 else
00151 {
00152 emean_[ieta0] = emean_[ieta0]+Original_Et;
00153 emean2[ieta0] = emean2[ieta0]+Original_Et*Original_Et;
00154 ntowers[ieta0]++;
00155 }
00156
00157 }
00158
00159 for(map<int,int>::const_iterator gt = geomtowers_.begin(); gt != geomtowers_.end(); gt++)
00160 {
00161
00162 int it = (*gt).first;
00163
00164 double e1 = (*(emean_.find(it))).second;
00165 double e2 = (*emean2.find(it)).second;
00166 int nt = (*gt).second - (*(ntowersWithJets_.find(it))).second;
00167
00168 LogDebug("PileUpSubtractor")<<" ieta : "<<it<<" number of towers : "<<nt<<" e1 : "<<e1<<" e2 : "<<e2<<"\n";
00169
00170 if(nt > 0) {
00171 emean_[it] = e1/nt;
00172 double eee = e2/nt - e1*e1/(nt*nt);
00173 if(eee<0.) eee = 0.;
00174 esigma_[it] = nSigmaPU_*sqrt(eee);
00175 }
00176 else
00177 {
00178 emean_[it] = 0.;
00179 esigma_[it] = 0.;
00180 }
00181 LogDebug("PileUpSubtractor")<<" ieta : "<<it<<" Pedestals : "<<emean_[it]<<" "<<esigma_[it]<<"\n";
00182 }
00183 }
00184
00185
00186
00187
00188
00189 double MultipleAlgoIterator::getEt(const reco::CandidatePtr & in) const {
00190 const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
00191 const GlobalPoint& pos=geo_->getPosition(ctc->id());
00192 double energy = ctc->emEnergy() + ctc->hadEnergy();
00193 double et = energy*sin(pos.theta());
00194 return et;
00195 }
00196
00197 double MultipleAlgoIterator::getEta(const reco::CandidatePtr & in) const {
00198 const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
00199 const GlobalPoint& pos=geo_->getPosition(ctc->id());
00200 double eta = pos.eta();
00201 return eta;
00202 }