CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_4/src/RecoHI/HiJetAlgos/src/ReflectedIterator.cc

Go to the documentation of this file.
00001 #include "RecoHI/HiJetAlgos/interface/ReflectedIterator.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 ReflectedIterator::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 ReflectedIterator::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 ReflectedIterator::subtractPedestal(vector<fastjet::PseudoJet> & coll)
00068 {
00069 
00070    LogDebug("PileUpSubtractor")<<"The subtractor subtracting pedestals...\n";
00071 
00072    int it = -100;
00073 
00074    vector<fastjet::PseudoJet> newcoll;
00075 
00076    for (vector<fastjet::PseudoJet>::iterator input_object = coll.begin (),
00077            fjInputsEnd = coll.end(); 
00078         input_object != fjInputsEnd; ++input_object) {
00079     
00080       reco::CandidatePtr const & itow =  (*inputs_)[ input_object->user_index() ];
00081     
00082       it = ieta( itow );
00083       iphi( itow );
00084     
00085       double Original_Et = itow->et();
00086       if(sumRecHits_){
00087          Original_Et = getEt(itow);
00088       }
00089 
00090       double etnew = Original_Et - (*(emean_.find(-it))).second - (*(esigma_.find(-it))).second;
00091       float mScale = etnew/input_object->Et(); 
00092       if(etnew < 0.) mScale = 0.;
00093     
00094       math::XYZTLorentzVectorD towP4(input_object->px()*mScale, input_object->py()*mScale,
00095                                      input_object->pz()*mScale, input_object->e()*mScale);
00096     
00097       int index = input_object->user_index();
00098       input_object->reset ( towP4.px(),
00099                             towP4.py(),
00100                             towP4.pz(),
00101                             towP4.energy() );
00102       input_object->set_user_index(index);
00103 
00104       if(etnew > 0. && dropZeroTowers_) newcoll.push_back(*input_object);
00105    }
00106 
00107    if(dropZeroTowers_) coll = newcoll;
00108 
00109 }
00110 
00111 void ReflectedIterator::calculatePedestal( vector<fastjet::PseudoJet> const & coll )
00112 {
00113    LogDebug("PileUpSubtractor")<<"The subtractor calculating pedestals...\n";
00114 
00115    map<int,double> emean2;
00116    map<int,int> ntowers;
00117     
00118    int ietaold = -10000;
00119    int ieta0 = -100;
00120    
00121    // Initial values for emean_, emean2, esigma_, ntowers
00122 
00123    for(int i = ietamin_; i < ietamax_+1; i++)
00124       {
00125          emean_[i] = 0.;
00126          emean2[i] = 0.;
00127          esigma_[i] = 0.;
00128          ntowers[i] = 0;
00129       }
00130     
00131    for (vector<fastjet::PseudoJet>::const_iterator input_object = coll.begin (),
00132            fjInputsEnd = coll.end();  
00133         input_object != fjInputsEnd; ++input_object) {
00134 
00135       const reco::CandidatePtr & originalTower=(*inputs_)[ input_object->user_index()];
00136       ieta0 = ieta( originalTower );
00137       double Original_Et = originalTower->et();
00138       if(sumRecHits_){
00139          Original_Et = getEt(originalTower);
00140       }
00141 
00142       if( ieta0-ietaold != 0 )
00143          {
00144             emean_[ieta0] = emean_[ieta0]+Original_Et;
00145             emean2[ieta0] = emean2[ieta0]+Original_Et*Original_Et;
00146             ntowers[ieta0] = 1;
00147             ietaold = ieta0;
00148          }
00149       else
00150          {
00151             emean_[ieta0] = emean_[ieta0]+Original_Et;
00152             emean2[ieta0] = emean2[ieta0]+Original_Et*Original_Et;
00153             ntowers[ieta0]++;
00154          }
00155 
00156    }
00157 
00158    for(map<int,int>::const_iterator gt = geomtowers_.begin(); gt != geomtowers_.end(); gt++)    
00159       {
00160 
00161          int it = (*gt).first;
00162        
00163          double e1 = (*(emean_.find(it))).second;
00164          double e2 = (*emean2.find(it)).second;
00165          int nt = (*gt).second - (*(ntowersWithJets_.find(it))).second;
00166 
00167          LogDebug("PileUpSubtractor")<<" ieta : "<<it<<" number of towers : "<<nt<<" e1 : "<<e1<<" e2 : "<<e2<<"\n";
00168         
00169          if(nt > 0) {
00170             emean_[it] = e1/nt;
00171             double eee = e2/nt - e1*e1/(nt*nt);    
00172             if(eee<0.) eee = 0.;
00173             esigma_[it] = nSigmaPU_*sqrt(eee);
00174          }
00175          else
00176             {
00177                emean_[it] = 0.;
00178                esigma_[it] = 0.;
00179             }
00180          LogDebug("PileUpSubtractor")<<" ieta : "<<it<<" Pedestals : "<<emean_[it]<<"  "<<esigma_[it]<<"\n";
00181       }
00182 }
00183 
00184 double ReflectedIterator::getEt(const reco::CandidatePtr & in) const {
00185    const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
00186    const GlobalPoint& pos=geo_->getPosition(ctc->id());
00187    double energy = ctc->emEnergy() + ctc->hadEnergy();
00188    double et = energy*sin(pos.theta());
00189    return et;
00190 }
00191 
00192 double ReflectedIterator::getEta(const reco::CandidatePtr & in) const {
00193    const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
00194    const GlobalPoint& pos=geo_->getPosition(ctc->id());
00195    double eta = pos.eta();
00196    return eta;
00197 }