CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_8_patch3/src/RecoHI/HiJetAlgos/src/MultipleAlgoIterator.cc

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