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
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 int ip = -100;
00111 vector<fastjet::PseudoJet> newcoll;
00112
00113 for (vector<fastjet::PseudoJet>::iterator input_object = coll.begin (),
00114 fjInputsEnd = coll.end();
00115 input_object != fjInputsEnd; ++input_object) {
00116
00117 reco::CandidatePtr const & itow = (*inputs_)[ input_object->user_index() ];
00118
00119 it = ieta( itow );
00120 ip = iphi( itow );
00121
00122 double Original_Et = itow->et();
00123 if(sumRecHits_){
00124 Original_Et = getEt(itow);
00125 }
00126
00127 double etnew = Original_Et - getPU(it,1,1);
00128 float mScale = etnew/input_object->Et();
00129 if(etnew < 0.) mScale = 0.;
00130
00131 math::XYZTLorentzVectorD towP4(input_object->px()*mScale, input_object->py()*mScale,
00132 input_object->pz()*mScale, input_object->e()*mScale);
00133
00134 int index = input_object->user_index();
00135 input_object->reset ( towP4.px(),
00136 towP4.py(),
00137 towP4.pz(),
00138 towP4.energy() );
00139 input_object->set_user_index(index);
00140 if(etnew > 0. && dropZeroTowers_) newcoll.push_back(*input_object);
00141 }
00142 if(dropZeroTowers_) coll = newcoll;
00143
00144 }
00145 }
00146
00147
00148 void ParametrizedSubtractor::calculateOrphanInput(vector<fastjet::PseudoJet> & orphanInput)
00149 {
00150 orphanInput = *fjInputs_;
00151 }
00152
00153
00154
00155 void ParametrizedSubtractor::offsetCorrectJets()
00156 {
00157
00158 LogDebug("PileUpSubtractor")<<"The subtractor correcting jets...\n";
00159 jetOffset_.clear();
00160
00161 using namespace reco;
00162
00163 (*fjInputs_) = fjOriginalInputs_;
00164 rescaleRMS(nSigmaPU_);
00165 subtractPedestal(*fjInputs_);
00166
00167 if(0){
00168 const fastjet::JetDefinition& def = fjClusterSeq_->jet_def();
00169 if ( !doAreaFastjet_ && !doRhoFastjet_) {
00170 fastjet::ClusterSequence newseq( *fjInputs_, def );
00171 (*fjClusterSeq_) = newseq;
00172 } else {
00173 fastjet::ClusterSequenceArea newseq( *fjInputs_, def , *fjActiveArea_ );
00174 (*fjClusterSeq_) = newseq;
00175 }
00176
00177 (*fjJets_) = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
00178 }
00179
00180 jetOffset_.reserve(fjJets_->size());
00181
00182 vector<fastjet::PseudoJet>::iterator pseudojetTMP = fjJets_->begin (),
00183 jetsEnd = fjJets_->end();
00184 for (; pseudojetTMP != jetsEnd; ++pseudojetTMP) {
00185
00186 int ijet = pseudojetTMP - fjJets_->begin();
00187 jetOffset_[ijet] = 0;
00188
00189 std::vector<fastjet::PseudoJet> towers =
00190 sorted_by_pt(fjClusterSeq_->constituents(*pseudojetTMP));
00191
00192 double newjetet = 0.;
00193 for(vector<fastjet::PseudoJet>::const_iterator ito = towers.begin(),
00194 towEnd = towers.end();
00195 ito != towEnd;
00196 ++ito)
00197 {
00198 const reco::CandidatePtr& originalTower = (*inputs_)[ito->user_index()];
00199 int it = ieta( originalTower );
00200 double Original_Et = originalTower->et();
00201
00202 if(sumRecHits_){
00203 Original_Et = getEt(originalTower);
00204 }
00205
00206 double etnew = Original_Et - getPU(it,1,1);
00207 if(etnew < 0.) etnew = 0;
00208 newjetet = newjetet + etnew;
00209 jetOffset_[ijet] += Original_Et - etnew;
00210 }
00211
00212 if(sumRecHits_){
00213 double mScale = newjetet/pseudojetTMP->Et();
00214 int cshist = pseudojetTMP->cluster_hist_index();
00215 pseudojetTMP->reset(pseudojetTMP->px()*mScale, pseudojetTMP->py()*mScale,
00216 pseudojetTMP->pz()*mScale, pseudojetTMP->e()*mScale);
00217 pseudojetTMP->set_cluster_hist_index(cshist);
00218 }
00219
00220 }
00221 }
00222
00223
00224 double ParametrizedSubtractor::getEt(const reco::CandidatePtr & in) const {
00225 const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
00226 const GlobalPoint& pos=geo_->getPosition(ctc->id());
00227 double energy = ctc->emEnergy() + ctc->hadEnergy();
00228
00229 if(0){
00230 energy = 0;
00231 const std::vector<DetId>& hitids = ctc->constituents();
00232 for(unsigned int i = 0; i< hitids.size(); ++i){
00233
00234 }
00235 }
00236
00237 double et = energy*sin(pos.theta());
00238 return et;
00239 }
00240
00241 double ParametrizedSubtractor::getEta(const reco::CandidatePtr & in) const {
00242 const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
00243 const GlobalPoint& pos=geo_->getPosition(ctc->id());
00244 double eta = pos.eta();
00245 return eta;
00246 }
00247
00248 double ParametrizedSubtractor::getMeanAtTower(const reco::CandidatePtr & in) const{
00249 int it = ieta(in);
00250 return getPU(it,1,0);
00251 }
00252
00253 double ParametrizedSubtractor::getSigmaAtTower(const reco::CandidatePtr & in) const {
00254 int it = ieta(in);
00255 return getPU(it,0,1);
00256 }
00257
00258 double ParametrizedSubtractor::getPileUpAtTower(const reco::CandidatePtr & in) const {
00259 int it = ieta(in);
00260 return getPU(it,1,1);
00261 }
00262
00263 double ParametrizedSubtractor::getPU(int ieta,bool addMean, bool addSigma) const {
00264
00265
00266
00267
00268 double em = hEtaMean[bin_]->GetBinContent(hEtaMean[bin_]->FindBin(ieta));
00269 double cm = fMean->Eval(centrality_);
00270
00271 double er = hEtaRMS[bin_]->GetBinContent(hEtaRMS[bin_]->FindBin(ieta));
00272 double cr = fRMS->Eval(centrality_);
00273
00274 if(interpolate_){
00275 double n = 0;
00276 int hbin = 40-bin_;
00277 double centerweight = (centrality_ - hC->GetBinCenter(hbin));
00278 double lowerweight = (centrality_ - hC->GetBinLowEdge(hbin));
00279 double upperweight = (centrality_ - hC->GetBinLowEdge(hbin+1));
00280
00281 em *= lowerweight*upperweight;
00282 er *= lowerweight*upperweight;
00283 n += lowerweight*upperweight;
00284
00285 if(bin_ > 0){
00286 em += upperweight*centerweight*hEtaMean[bin_]->GetBinContent(hEtaMean[bin_-1]->FindBin(ieta));
00287 er += upperweight*centerweight*hEtaRMS[bin_]->GetBinContent(hEtaRMS[bin_-1]->FindBin(ieta));
00288 n += upperweight*centerweight;
00289 }
00290
00291 if(bin_ < 39){
00292 em += lowerweight*centerweight*hEtaMean[bin_]->GetBinContent(hEtaMean[bin_+1]->FindBin(ieta));
00293 er += lowerweight*centerweight*hEtaRMS[bin_]->GetBinContent(hEtaRMS[bin_+1]->FindBin(ieta));
00294 n += lowerweight*centerweight;
00295 }
00296 em /= n;
00297 er /= n;
00298 }
00299
00300
00301 return addMean*em*cm + addSigma*nSigmaPU_*er*cr;
00302 }
00303
00304