CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/RecoHI/HiJetAlgos/src/ParametrizedSubtractor.cc

Go to the documentation of this file.
00001 #include "RecoHI/HiJetAlgos/interface/ParametrizedSubtractor.h"
00002 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00003 #include "DataFormats/Candidate/interface/Candidate.h"
00004 #include "FWCore/ParameterSet/interface/FileInPath.h"
00005 #include "FWCore/Framework/interface/ESHandle.h"
00006 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
00007 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00008 
00009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00010 #include "string"
00011 #include "iostream"
00012 using namespace std;
00013 
00014 void ParametrizedSubtractor::rescaleRMS(double s){
00015    for ( std::map<int, double>::iterator iter = esigma_.begin();
00016          iter != esigma_.end(); ++iter ){
00017       iter->second = s*(iter->second);
00018    }
00019 }
00020 
00021 
00022 ParametrizedSubtractor::ParametrizedSubtractor(const edm::ParameterSet& iConfig) : 
00023    PileUpSubtractor(iConfig),
00024    dropZeroTowers_(iConfig.getUntrackedParameter<bool>("dropZeroTowers",true)),
00025    cbins_(0)
00026 {
00027    centTag_ = iConfig.getUntrackedParameter<edm::InputTag>("centTag",edm::InputTag("hiCentrality","","RECO"));
00028 
00029    interpolate_ = iConfig.getParameter<bool>("interpolate");
00030    sumRecHits_ = iConfig.getParameter<bool>("sumRecHits");
00031 
00032    std::string ifname = "RecoHI/HiJetAlgos/data/PU_DATA.root";
00033    TFile* inf = new TFile(edm::FileInPath(ifname).fullPath().data());
00034    fPU = (TF1*)inf->Get("fPU");
00035    fMean = (TF1*)inf->Get("fMean");
00036    fRMS = (TF1*)inf->Get("fRMS");
00037    hC = (TH1D*)inf->Get("hC");
00038 
00039    for(int i = 0; i < 40; ++i){
00040       hEta.push_back((TH1D*)inf->Get(Form("hEta_%d",i)));
00041       hEtaMean.push_back((TH1D*)inf->Get(Form("hEtaMean_%d",i)));
00042       hEtaRMS.push_back((TH1D*)inf->Get(Form("hEtaRMS_%d",i)));
00043    }
00044 
00045 }
00046 
00047 
00048 void ParametrizedSubtractor::setupGeometryMap(edm::Event& iEvent,const edm::EventSetup& iSetup){
00049    LogDebug("PileUpSubtractor")<<"The subtractor setting up geometry...\n";
00050 
00051    //   if(!cbins_) getCentralityBinsFromDB(iSetup);
00052 
00053    edm::Handle<reco::Centrality> cent;
00054    iEvent.getByLabel(centTag_,cent);
00055    
00056    centrality_ = cent->EtHFhitSum();
00057    bin_ = 40-hC->FindBin(centrality_);
00058    if(bin_ > 39) bin_ = 39;
00059    if(bin_ < 0) bin_ = 0;
00060 
00061    if(!geo_) {
00062       edm::ESHandle<CaloGeometry> pG;
00063       iSetup.get<CaloGeometryRecord>().get(pG);
00064       geo_ = pG.product();
00065       std::vector<DetId> alldid =  geo_->getValidDetIds();
00066       
00067       int ietaold = -10000;
00068       ietamax_ = -10000;
00069       ietamin_ = 10000;
00070       for(std::vector<DetId>::const_iterator did=alldid.begin(); did != alldid.end(); did++){
00071          if( (*did).det() == DetId::Hcal ){
00072             HcalDetId hid = HcalDetId(*did);
00073             if( (hid).depth() == 1 ) {
00074                allgeomid_.push_back(*did);
00075 
00076                if((hid).ieta() != ietaold){
00077                   ietaold = (hid).ieta();
00078                   geomtowers_[(hid).ieta()] = 1;
00079                   if((hid).ieta() > ietamax_) ietamax_ = (hid).ieta();
00080                   if((hid).ieta() < ietamin_) ietamin_ = (hid).ieta();
00081                }
00082                else{
00083                   geomtowers_[(hid).ieta()]++;
00084                }
00085             }
00086          }
00087       }
00088    }
00089 
00090    for (int i = ietamin_; i<ietamax_+1; i++) {
00091       emean_[i] = 0.;
00092       esigma_[i] = 0.;
00093       ntowersWithJets_[i] = 0;
00094    }
00095 }
00096 
00097 
00098 void ParametrizedSubtractor::calculatePedestal( vector<fastjet::PseudoJet> const & coll ){
00099    return;
00100 }
00101 
00102 void ParametrizedSubtractor::subtractPedestal(vector<fastjet::PseudoJet> & coll)
00103 {
00104    if(0){
00105       return;
00106    }else{
00107       LogDebug("PileUpSubtractor")<<"The subtractor subtracting pedestals...\n";
00108 
00109    int it = -100;
00110    vector<fastjet::PseudoJet> newcoll;
00111         
00112    for (vector<fastjet::PseudoJet>::iterator input_object = coll.begin (),
00113            fjInputsEnd = coll.end(); 
00114         input_object != fjInputsEnd; ++input_object) {
00115     
00116       reco::CandidatePtr const & itow =  (*inputs_)[ input_object->user_index() ];
00117     
00118       it = ieta( itow );
00119       iphi( itow );
00120 
00121       double Original_Et = itow->et();
00122       if(sumRecHits_){
00123          Original_Et = getEt(itow);
00124       }
00125 
00126       double etnew = Original_Et - getPU(it,1,1);
00127       float mScale = etnew/input_object->Et(); 
00128       if(etnew < 0.) mScale = 0.;
00129     
00130       math::XYZTLorentzVectorD towP4(input_object->px()*mScale, input_object->py()*mScale,
00131                                      input_object->pz()*mScale, input_object->e()*mScale);
00132     
00133       int index = input_object->user_index();
00134       input_object->reset ( towP4.px(),
00135                             towP4.py(),
00136                             towP4.pz(),
00137                             towP4.energy() );
00138       input_object->set_user_index(index);
00139       if(etnew > 0. && dropZeroTowers_) newcoll.push_back(*input_object);
00140    }
00141    if(dropZeroTowers_) coll = newcoll;
00142 
00143    }
00144 }
00145 
00146 
00147 void ParametrizedSubtractor::calculateOrphanInput(vector<fastjet::PseudoJet> & orphanInput) 
00148 {
00149    orphanInput = *fjInputs_;
00150 }
00151 
00152 
00153 
00154 void ParametrizedSubtractor::offsetCorrectJets()
00155 {
00156 
00157   LogDebug("PileUpSubtractor")<<"The subtractor correcting jets...\n";
00158   jetOffset_.clear();
00159 
00160   using namespace reco;
00161   
00162   (*fjInputs_) = fjOriginalInputs_;
00163   rescaleRMS(nSigmaPU_);
00164   subtractPedestal(*fjInputs_);
00165 
00166   if(0){
00167      const fastjet::JetDefinition& def = fjClusterSeq_->jet_def();
00168      if ( !doAreaFastjet_ && !doRhoFastjet_) {
00169         fastjet::ClusterSequence newseq( *fjInputs_, def );
00170         (*fjClusterSeq_) = newseq;
00171      } else {
00172         fastjet::ClusterSequenceArea newseq( *fjInputs_, def , *fjActiveArea_ );
00173         (*fjClusterSeq_) = newseq;
00174      }
00175      
00176      (*fjJets_) = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
00177   }
00178   
00179   jetOffset_.reserve(fjJets_->size());
00180   
00181   vector<fastjet::PseudoJet>::iterator pseudojetTMP = fjJets_->begin (),
00182     jetsEnd = fjJets_->end();
00183   for (; pseudojetTMP != jetsEnd; ++pseudojetTMP) {
00184     
00185     int ijet = pseudojetTMP - fjJets_->begin();
00186     jetOffset_[ijet] = 0;
00187     
00188     std::vector<fastjet::PseudoJet> towers =
00189       sorted_by_pt(fjClusterSeq_->constituents(*pseudojetTMP));
00190     
00191     double newjetet = 0.;
00192     for(vector<fastjet::PseudoJet>::const_iterator ito = towers.begin(),
00193           towEnd = towers.end();
00194         ito != towEnd;
00195         ++ito)
00196       {
00197         const reco::CandidatePtr& originalTower = (*inputs_)[ito->user_index()];
00198         int it = ieta( originalTower );
00199         double Original_Et = originalTower->et();
00200       
00201         if(sumRecHits_){
00202            Original_Et = getEt(originalTower);
00203         }
00204 
00205         double etnew = Original_Et - getPU(it,1,1);
00206         if(etnew < 0.) etnew = 0;
00207         newjetet = newjetet + etnew;
00208         jetOffset_[ijet] += Original_Et - etnew;
00209       }
00210 
00211     if(sumRecHits_){       
00212        double mScale = newjetet/pseudojetTMP->Et();
00213        int cshist = pseudojetTMP->cluster_hist_index();
00214        pseudojetTMP->reset(pseudojetTMP->px()*mScale, pseudojetTMP->py()*mScale,
00215                            pseudojetTMP->pz()*mScale, pseudojetTMP->e()*mScale);
00216        pseudojetTMP->set_cluster_hist_index(cshist);
00217     }
00218 
00219   }
00220 }
00221 
00222 
00223 double ParametrizedSubtractor::getEt(const reco::CandidatePtr & in) const {
00224    const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
00225    const GlobalPoint& pos=geo_->getPosition(ctc->id());
00226    double energy = ctc->emEnergy() + ctc->hadEnergy();
00227 
00228    if(0){
00229       energy = 0;
00230       const std::vector<DetId>& hitids = ctc->constituents();
00231       for(unsigned int i = 0; i< hitids.size(); ++i){
00232 
00233       }
00234    }
00235 
00236    double et = energy*sin(pos.theta());
00237    return et;
00238 }
00239 
00240 double ParametrizedSubtractor::getEta(const reco::CandidatePtr & in) const {
00241    const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
00242    const GlobalPoint& pos=geo_->getPosition(ctc->id());
00243    double eta = pos.eta();
00244    return eta;
00245 }
00246 
00247 double ParametrizedSubtractor::getMeanAtTower(const reco::CandidatePtr & in) const{
00248    int it = ieta(in);
00249    return getPU(it,1,0);
00250 }
00251 
00252 double ParametrizedSubtractor::getSigmaAtTower(const reco::CandidatePtr & in) const {
00253    int it = ieta(in);
00254    return getPU(it,0,1);
00255 }
00256 
00257 double ParametrizedSubtractor::getPileUpAtTower(const reco::CandidatePtr & in) const {
00258    int it = ieta(in);
00259    return getPU(it,1,1);
00260 }
00261 
00262 double ParametrizedSubtractor::getPU(int ieta,bool addMean, bool addSigma) const {   
00263 
00264   //double e = hEta[bin_]->GetBinContent(hEta[bin_]->FindBin(ieta));
00265   //double c = fPU->Eval(centrality_);
00266 
00267    double em = hEtaMean[bin_]->GetBinContent(hEtaMean[bin_]->FindBin(ieta));
00268    double cm = fMean->Eval(centrality_);
00269 
00270    double er = hEtaRMS[bin_]->GetBinContent(hEtaRMS[bin_]->FindBin(ieta));
00271    double cr = fRMS->Eval(centrality_);
00272 
00273    if(interpolate_){
00274       double n = 0;
00275       int hbin = 40-bin_;
00276       double centerweight =  (centrality_ - hC->GetBinCenter(hbin));
00277       double lowerweight = (centrality_ - hC->GetBinLowEdge(hbin));
00278       double upperweight = (centrality_ - hC->GetBinLowEdge(hbin+1));
00279 
00280       em *= lowerweight*upperweight;
00281       er *= lowerweight*upperweight;
00282       n += lowerweight*upperweight;
00283 
00284       if(bin_ > 0){
00285          em += upperweight*centerweight*hEtaMean[bin_]->GetBinContent(hEtaMean[bin_-1]->FindBin(ieta));
00286          er += upperweight*centerweight*hEtaRMS[bin_]->GetBinContent(hEtaRMS[bin_-1]->FindBin(ieta));
00287          n += upperweight*centerweight;
00288       }
00289 
00290       if(bin_ < 39){
00291          em += lowerweight*centerweight*hEtaMean[bin_]->GetBinContent(hEtaMean[bin_+1]->FindBin(ieta));
00292          er += lowerweight*centerweight*hEtaRMS[bin_]->GetBinContent(hEtaRMS[bin_+1]->FindBin(ieta));
00293          n += lowerweight*centerweight;
00294       }
00295       em /= n;
00296       er /= n;
00297    }
00298 
00299    //   return e*c;
00300    return addMean*em*cm + addSigma*nSigmaPU_*er*cr;
00301 }
00302 
00303