00001
00002 #include <memory>
00003 #include <string>
00004
00005
00006 #include "FWCore/Framework/interface/Frameworkfwd.h"
00007 #include "FWCore/Framework/interface/EDAnalyzer.h"
00008
00009 #include "FWCore/Framework/interface/Event.h"
00010 #include "FWCore/Framework/interface/MakerMacros.h"
00011
00012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00013
00014 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00015 #include "FWCore/ServiceRegistry/interface/Service.h"
00016
00017 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00018
00019 #include "AnalysisDataFormats/TopObjects/interface/TtGenEvent.h"
00020
00021 #include "DataFormats/PatCandidates/interface/Muon.h"
00022 #include "DataFormats/PatCandidates/interface/Jet.h"
00023 #include "DataFormats/PatCandidates/interface/MET.h"
00024 #include "DataFormats/PatCandidates/interface/Tau.h"
00025 #include "DataFormats/PatCandidates/interface/Electron.h"
00026
00027 #include "TDirectory.h"
00028 #include "TH1F.h"
00029 #include "TF1.h"
00030 #include "TTree.h"
00031
00032 #include <Math/VectorUtil.h>
00033
00034
00035
00036
00037
00038 class ResolutionCreator : public edm::EDAnalyzer {
00039 public:
00040 explicit ResolutionCreator(const edm::ParameterSet&);
00041 ~ResolutionCreator();
00042
00043 private:
00044 virtual void beginJob() ;
00045 virtual void analyze(const edm::Event&, const edm::EventSetup&);
00046 virtual void endJob() ;
00047
00048
00049 std::string objectType_, labelName_;
00050 std::vector<double> etabinVals_, pTbinVals_;
00051 double minDR_;
00052 int matchingAlgo_;
00053 bool useMaxDist_;
00054 bool useDeltaR_;
00055 double maxDist_;
00056 int ptnrbins, etanrbins;
00057 int nrFilled;
00058
00059
00060 TF1 *fResPtEtaBin[10][20][20];
00061 TF1 *fResEtaBin[10][20];
00062 TH1F *hResPtEtaBin[10][20][20];
00063 TH1F *hResEtaBin[10][20];
00064 TTree* tResVar;
00065 };
00066
00067
00068
00069
00070
00071 ResolutionCreator::ResolutionCreator(const edm::ParameterSet& iConfig)
00072 {
00073
00074 objectType_ = iConfig.getParameter< std::string > ("object");
00075 labelName_ = iConfig.getParameter< std::string > ("label");
00076 if(objectType_ != "met"){
00077 etabinVals_ = iConfig.getParameter< std::vector<double> > ("etabinValues");
00078 }
00079 pTbinVals_ = iConfig.getParameter< std::vector<double> > ("pTbinValues");
00080 minDR_ = iConfig.getParameter< double > ("minMatchingDR");
00081
00082 nrFilled = 0;
00083
00084 }
00085
00086
00087 ResolutionCreator::~ResolutionCreator()
00088 {
00089 }
00090
00091
00092
00093
00094
00095
00096
00097 void
00098 ResolutionCreator::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00099 {
00100
00101
00102 std::vector<reco::Particle *> p4gen, p4rec;
00103
00104 edm::Handle<TtGenEvent> genEvt;
00105 iEvent.getByLabel ("genEvt",genEvt);
00106
00107 if(genEvt->particles().size()<10) return;
00108
00109 if(objectType_ == "electron"){
00110 edm::Handle<std::vector<pat::Electron> > electrons;
00111 iEvent.getByLabel(labelName_,electrons);
00112 for(size_t e=0; e<electrons->size(); e++) {
00113 for(size_t p=0; p<genEvt->particles().size(); p++){
00114 if( (std::abs(genEvt->particles()[p].pdgId()) == 11) && (ROOT::Math::VectorUtil::DeltaR(genEvt->particles()[p].p4(), (*electrons)[e].p4()) < minDR_) ) {
00115
00116
00117 }
00118 }
00119 }
00120 }
00121 else if(objectType_ == "muon"){
00122 edm::Handle<std::vector<pat::Muon> > muons;
00123 iEvent.getByLabel(labelName_,muons);
00124 for(size_t m=0; m<muons->size(); m++) {
00125 for(size_t p=0; p<genEvt->particles().size(); p++){
00126 if( (std::abs(genEvt->particles()[p].pdgId()) == 13) && (ROOT::Math::VectorUtil::DeltaR(genEvt->particles()[p].p4(), (*muons)[m].p4()) < minDR_) ) {
00127
00128
00129 }
00130 }
00131 }
00132 }
00133 else if(objectType_ == "lJets" ){
00134 edm::Handle<std::vector<pat::Jet> > jets;
00135 iEvent.getByLabel(labelName_,jets);
00136 if(jets->size()>=4) {
00137 for(unsigned int j = 0; j<4; j++){
00138 for(size_t p=0; p<genEvt->particles().size(); p++){
00139 if( (std::abs(genEvt->particles()[p].pdgId()) < 5) && (ROOT::Math::VectorUtil::DeltaR(genEvt->particles()[p].p4(), (*jets)[j].p4())< minDR_) ){
00140
00141
00142 }
00143 }
00144 }
00145 }
00146 }
00147 else if(objectType_ == "bJets" ){
00148 edm::Handle<std::vector<pat::Jet> > jets;
00149 iEvent.getByLabel(labelName_,jets);
00150 if(jets->size()>=4) {
00151 for(unsigned int j = 0; j<4; j++){
00152 for(size_t p=0; p<genEvt->particles().size(); p++){
00153 if( (std::abs(genEvt->particles()[p].pdgId()) == 5) && (ROOT::Math::VectorUtil::DeltaR(genEvt->particles()[p].p4(), (*jets)[j].p4())< minDR_) ) {
00154
00155
00156 }
00157 }
00158 }
00159 }
00160 }
00161 else if(objectType_ == "met"){
00162 edm::Handle<std::vector<pat::MET> > mets;
00163 iEvent.getByLabel(labelName_,mets);
00164 if(mets->size()>=1) {
00165 if( genEvt->isSemiLeptonic() && genEvt->singleNeutrino() != 0 && ROOT::Math::VectorUtil::DeltaR(genEvt->singleNeutrino()->p4(), (*mets)[0].p4()) < minDR_) {
00166
00167
00168 }
00169 }
00170 }
00171 else if(objectType_ == "tau"){
00172 edm::Handle<std::vector<pat::Tau> > taus;
00173 iEvent.getByLabel(labelName_,taus);
00174 for(std::vector<pat::Tau>::const_iterator tau = taus->begin(); tau != taus->end(); ++tau) {
00175
00176 reco::GenParticle genLepton = *(tau->genLepton());
00177 if( std::abs(genLepton.pdgId())==15 && genLepton.status()==2 &&
00178 genLepton.numberOfMothers()>0 &&
00179 std::abs(genLepton.mother(0)->pdgId())==15 &&
00180 genLepton.mother(0)->numberOfMothers()>0 &&
00181 std::abs(genLepton.mother(0)->mother(0)->pdgId())==24 &&
00182 ROOT::Math::VectorUtil::DeltaR(genLepton.p4(), tau->p4()) < minDR_ ) {
00183 }
00184
00185
00186 }
00187 }
00188
00189 for(unsigned m=0; m<p4gen.size(); m++){
00190 double Egen = p4gen[m]->energy();
00191 double Thetagen = p4gen[m]->theta();
00192 double Phigen = p4gen[m]->phi();
00193 double Etgen = p4gen[m]->et();
00194 double Etagen = p4gen[m]->eta();
00195 double Ecal = p4rec[m]->energy();
00196 double Thetacal = p4rec[m]->theta();
00197 double Phical = p4rec[m]->phi();
00198 double Etcal = p4rec[m]->et();
00199 double Etacal = p4rec[m]->eta();
00200 double phidiff = Phical- Phigen;
00201 if(phidiff>3.14159) phidiff = 2.*3.14159 - phidiff;
00202 if(phidiff<-3.14159) phidiff = -phidiff - 2.*3.14159;
00203
00204
00205 int etabin = 0;
00206 if(etanrbins > 1){
00207 for(unsigned int b=0; b<etabinVals_.size()-1; b++) {
00208 if(fabs(Etacal) > etabinVals_[b]) etabin = b;
00209 }
00210 }
00211
00212 int ptbin = 0;
00213 for(unsigned int b=0; b<pTbinVals_.size()-1; b++) {
00214 if(p4rec[m]->pt() > pTbinVals_[b]) ptbin = b;
00215 }
00216
00217
00218
00219
00220
00221
00222
00223
00224
00225
00226
00227
00228 ROOT::Math::SVector<double,3> pcalvec(p4rec[m]->px(),p4rec[m]->py(),p4rec[m]->pz());
00229 ROOT::Math::SVector<double,3> pgenvec(p4gen[m]->px(),p4gen[m]->py(),p4gen[m]->pz());
00230
00231 ROOT::Math::SVector<double,3> u_z(0,0,1);
00232 ROOT::Math::SVector<double,3> u_1 = ROOT::Math::Unit(pcalvec);
00233 ROOT::Math::SVector<double,3> u_3 = ROOT::Math::Cross(u_z,u_1)/ROOT::Math::Mag(ROOT::Math::Cross(u_z,u_1));
00234 ROOT::Math::SVector<double,3> u_2 = ROOT::Math::Cross(u_1,u_3)/ROOT::Math::Mag(ROOT::Math::Cross(u_1,u_3));
00235 double acal = 1.;
00236 double bcal = 0.;
00237 double ccal = 0.;
00238 double dcal = 1.;
00239 double agen = ROOT::Math::Dot(pgenvec,u_1)/ROOT::Math::Mag(pcalvec);
00240 double bgen = ROOT::Math::Dot(pgenvec,u_2);
00241 double cgen = ROOT::Math::Dot(pgenvec,u_3);
00242 double dgen = Egen/Ecal;
00243
00244
00245 ++nrFilled;
00246 hResPtEtaBin[0][etabin][ptbin] -> Fill(acal-agen);
00247 hResPtEtaBin[1][etabin][ptbin] -> Fill(bcal-bgen);
00248 hResPtEtaBin[2][etabin][ptbin] -> Fill(ccal-cgen);
00249 hResPtEtaBin[3][etabin][ptbin] -> Fill(dcal-dgen);
00250 hResPtEtaBin[4][etabin][ptbin] -> Fill(Thetacal-Thetagen);
00251 hResPtEtaBin[5][etabin][ptbin] -> Fill(phidiff);
00252 hResPtEtaBin[6][etabin][ptbin] -> Fill(Etcal-Etgen);
00253 hResPtEtaBin[7][etabin][ptbin] -> Fill(Etacal-Etagen);
00254
00255 delete p4gen[m];
00256 delete p4rec[m];
00257 }
00258
00259 }
00260
00261
00262
00263 void
00264 ResolutionCreator::beginJob()
00265 {
00266 edm::Service<TFileService> fs;
00267 if (!fs) throw edm::Exception(edm::errors::Configuration, "TFileService missing from configuration!");
00268
00269
00270 TString resObsName[8] = {"_ares","_bres","_cres","_dres","_thres","_phres","_etres","_etares"};
00271 int resObsNrBins = 120;
00272 if( (objectType_ == "muon") || (objectType_ == "electron") ) resObsNrBins = 80;
00273 std::vector<double> resObsMin, resObsMax;
00274 if(objectType_ == "electron"){
00275 resObsMin.push_back(-0.15); resObsMin.push_back(-0.2); resObsMin.push_back(-0.1); resObsMin.push_back(-0.15); resObsMin.push_back(-0.0012); resObsMin.push_back(-0.009); resObsMin.push_back(-16); resObsMin.push_back(-0.0012);
00276 resObsMax.push_back( 0.15); resObsMax.push_back( 0.2); resObsMax.push_back( 0.1); resObsMax.push_back( 0.15); resObsMax.push_back( 0.0012); resObsMax.push_back( 0.009); resObsMax.push_back( 16); resObsMax.push_back( 0.0012);
00277 } else if(objectType_ == "muon"){
00278 resObsMin.push_back(-0.15); resObsMin.push_back(-0.1); resObsMin.push_back(-0.05); resObsMin.push_back(-0.15); resObsMin.push_back(-0.004); resObsMin.push_back(-0.003); resObsMin.push_back(-8); resObsMin.push_back(-0.004);
00279 resObsMax.push_back( 0.15); resObsMax.push_back( 0.1); resObsMax.push_back( 0.05); resObsMax.push_back( 0.15); resObsMax.push_back( 0.004); resObsMax.push_back( 0.003); resObsMax.push_back( 8); resObsMax.push_back( 0.004);
00280 } else if(objectType_ == "tau"){
00281 resObsMin.push_back(-1.); resObsMin.push_back(-10.); resObsMin.push_back(-10); resObsMin.push_back(-1.); resObsMin.push_back(-0.1); resObsMin.push_back(-0.1); resObsMin.push_back(-80); resObsMin.push_back(-0.1);
00282 resObsMax.push_back( 1.); resObsMax.push_back( 10.); resObsMax.push_back( 10); resObsMax.push_back( 1.); resObsMax.push_back( 0.1); resObsMax.push_back( 0.1); resObsMax.push_back( 50); resObsMax.push_back( 0.1);
00283 } else if(objectType_ == "lJets" || objectType_ == "bJets"){
00284 resObsMin.push_back(-1.); resObsMin.push_back(-10.); resObsMin.push_back(-10.); resObsMin.push_back(-1.); resObsMin.push_back(-0.4); resObsMin.push_back(-0.6); resObsMin.push_back( -80); resObsMin.push_back(-0.6);
00285 resObsMax.push_back( 1.); resObsMax.push_back( 10.); resObsMax.push_back( 10.); resObsMax.push_back( 1.); resObsMax.push_back( 0.4); resObsMax.push_back( 0.6); resObsMax.push_back( 80); resObsMax.push_back( 0.6);
00286 } else{
00287 resObsMin.push_back(-2.); resObsMin.push_back(-150.); resObsMin.push_back(-150.); resObsMin.push_back(-2.); resObsMin.push_back(-6); resObsMin.push_back(-6); resObsMin.push_back( -180); resObsMin.push_back(-6);
00288 resObsMax.push_back( 3.); resObsMax.push_back( 150.); resObsMax.push_back( 150.); resObsMax.push_back( 3.); resObsMax.push_back( 6); resObsMax.push_back( 6); resObsMax.push_back( 180); resObsMax.push_back( 6);
00289 }
00290
00291 const char* resObsVsPtFit[8] = { "[0]+[1]*exp(-[2]*x)",
00292 "[0]+[1]*exp(-[2]*x)",
00293 "[0]+[1]*exp(-[2]*x)",
00294 "[0]+[1]*exp(-[2]*x)",
00295 "[0]+[1]*exp(-[2]*x)",
00296 "[0]+[1]*exp(-[2]*x)",
00297 "pol1",
00298 "[0]+[1]*exp(-[2]*x)"
00299 };
00300
00301 ptnrbins = pTbinVals_.size()-1;
00302 double *ptbins = new double[pTbinVals_.size()];
00303 for(unsigned int b=0; b<pTbinVals_.size(); b++) ptbins[b] = pTbinVals_[b];
00304 double *etabins;
00305 if(objectType_ != "met"){
00306 etanrbins = etabinVals_.size()-1;
00307 etabins = new double[etabinVals_.size()];
00308 for(unsigned int b=0; b<etabinVals_.size(); b++) etabins[b] = etabinVals_[b];
00309 }else{
00310 etanrbins = 1;
00311 etabins = new double[2];
00312 etabins[0] = 0; etabins[1] = 5.;
00313 }
00314
00315
00316
00317 for(Int_t ro=0; ro<8; ro++) {
00318 for(Int_t etab=0; etab<etanrbins; etab++) {
00319 for(Int_t ptb=0; ptb<ptnrbins; ptb++) {
00320 TString obsName = objectType_; obsName += resObsName[ro]; obsName += "_etabin"; obsName += etab; obsName += "_ptbin";
00321 obsName += ptb;
00322 hResPtEtaBin[ro][etab][ptb] = fs->make<TH1F>(obsName,obsName,resObsNrBins,resObsMin[ro],resObsMax[ro]);
00323 fResPtEtaBin[ro][etab][ptb] = fs->make<TF1>("F_"+obsName,"gaus");
00324 }
00325 TString obsName2 = objectType_; obsName2 += resObsName[ro]; obsName2 += "_etabin"; obsName2 += etab;
00326 hResEtaBin[ro][etab] = fs->make<TH1F>(obsName2,obsName2,ptnrbins,ptbins);
00327 fResEtaBin[ro][etab] = fs->make<TF1>("F_"+obsName2,resObsVsPtFit[ro],pTbinVals_[0],pTbinVals_[pTbinVals_.size()-1]);
00328 }
00329 }
00330 tResVar = fs->make< TTree >("tResVar","Resolution tree");
00331
00332 delete [] etabins;
00333 delete [] ptbins;
00334
00335 }
00336
00337
00338 void
00339 ResolutionCreator::endJob() {
00340 TString resObsName2[8] = {"a","b","c","d","theta","phi","et","eta"};
00341 Int_t ro=0;
00342 Double_t pt=0.;
00343 Double_t eta=0.;
00344 Double_t value,error;
00345
00346 tResVar->Branch("Pt",&pt,"Pt/D");
00347 tResVar->Branch("Eta",&eta,"Eta/D");
00348 tResVar->Branch("ro",&ro,"ro/I");
00349 tResVar->Branch("value",&value,"value/D");
00350 tResVar->Branch("error",&error,"error/D");
00351
00352 for(ro=0; ro<8; ro++) {
00353 for(int etab=0; etab<etanrbins; etab++) {
00354
00355 eta = etanrbins > 1 ? (etabinVals_[etab]+etabinVals_[etab+1])/2. : 2.5 ;
00356 for(int ptb=0; ptb<ptnrbins; ptb++) {
00357
00358 pt = (pTbinVals_[ptb]+pTbinVals_[ptb+1])/2.;
00359 double maxcontent = 0.;
00360 int maxbin = 0;
00361 for(int nb=1; nb<hResPtEtaBin[ro][etab][ptb]->GetNbinsX(); nb ++){
00362 if (hResPtEtaBin[ro][etab][ptb]->GetBinContent(nb)>maxcontent) {
00363 maxcontent = hResPtEtaBin[ro][etab][ptb]->GetBinContent(nb);
00364 maxbin = nb;
00365 }
00366 }
00367 int range = (int)(hResPtEtaBin[ro][etab][ptb]->GetNbinsX()/6);
00368 fResPtEtaBin[ro][etab][ptb] -> SetRange(hResPtEtaBin[ro][etab][ptb]->GetBinCenter(maxbin-range),
00369 hResPtEtaBin[ro][etab][ptb]->GetBinCenter(maxbin+range));
00370 fResPtEtaBin[ro][etab][ptb] -> SetParameters(hResPtEtaBin[ro][etab][ptb] -> GetMaximum(),
00371 hResPtEtaBin[ro][etab][ptb] -> GetMean(),
00372 hResPtEtaBin[ro][etab][ptb] -> GetRMS());
00373 hResPtEtaBin[ro][etab][ptb] -> Fit(fResPtEtaBin[ro][etab][ptb]->GetName(),"RQ");
00374 hResEtaBin[ro][etab] -> SetBinContent(ptb+1,fResPtEtaBin[ro][etab][ptb]->GetParameter(2));
00375 hResEtaBin[ro][etab] -> SetBinError(ptb+1,fResPtEtaBin[ro][etab][ptb]->GetParError(2));
00376
00377 value = fResPtEtaBin[ro][etab][ptb]->GetParameter(2);
00378 error = fResPtEtaBin[ro][etab][ptb]->GetParError(2);
00379 tResVar->Fill();
00380 }
00381
00382
00383 pt = 0.;
00384 value = ((pTbinVals_[0]+pTbinVals_[1])/2.)*(fResPtEtaBin[ro][etab][0]->GetParameter(2)-fResPtEtaBin[ro][etab][1]->GetParameter(2))/((pTbinVals_[2]-pTbinVals_[0])/2.)+fResPtEtaBin[ro][etab][0]->GetParameter(2);
00385 error = fResPtEtaBin[ro][etab][0]->GetParError(2)+fResPtEtaBin[ro][etab][1]->GetParError(2);
00386 tResVar->Fill();
00387
00388 hResEtaBin[ro][etab] -> Fit(fResEtaBin[ro][etab]->GetName(),"RQ");
00389 }
00390 }
00391 if(objectType_ == "lJets" && nrFilled == 0) edm::LogProblem ("SummaryError") << "No plots filled for light jets \n";
00392 if(objectType_ == "bJets" && nrFilled == 0) edm::LogProblem ("SummaryError") << "No plots filled for bjets \n";
00393 if(objectType_ == "muon" && nrFilled == 0) edm::LogProblem ("SummaryError") << "No plots filled for muons \n";
00394 if(objectType_ == "electron" && nrFilled == 0) edm::LogProblem ("SummaryError") << "No plots filled for electrons \n";
00395 if(objectType_ == "tau" && nrFilled == 0) edm::LogProblem ("SummaryError") << "No plots filled for taus \n";
00396 if(objectType_ == "met" && nrFilled == 0) edm::LogProblem ("SummaryError") << "No plots filled for met \n";
00397
00398 edm::LogVerbatim ("MainResults") << " \n\n";
00399 edm::LogVerbatim ("MainResults") << " ----------------------------------------------";
00400 edm::LogVerbatim ("MainResults") << " ----------------------------------------------";
00401 edm::LogVerbatim ("MainResults") << " Resolutions on "<< objectType_ << " with nrfilled: "<<nrFilled;
00402 edm::LogVerbatim ("MainResults") << " ----------------------------------------------";
00403 edm::LogVerbatim ("MainResults") << " ----------------------------------------------";
00404 if(nrFilled != 0 && objectType_ != "met") {
00405 for(ro=0; ro<8; ro++) {
00406 edm::LogVerbatim ("MainResults") << "-------------------- ";
00407 edm::LogVerbatim ("MainResults") << "\n Resolutions on " << resObsName2[ro] << "\n";
00408 edm::LogVerbatim ("MainResults") << "-------------------- ";
00409 for(int etab=0; etab<etanrbins; etab++) {
00410 if(nrFilled != 0 && ro != 6) {
00411 if(etab == 0){
00412 edm::LogVerbatim ("MainResults") << "if(fabs(eta)<"<<etabinVals_[etab+1] <<") res = " <<
00413 fResEtaBin[ro][etab]->GetParameter(0) << "+" << fResEtaBin[ro][etab]->GetParameter(1)
00414 << "*exp(-(" <<fResEtaBin[ro][etab]->GetParameter(2) << "*pt));";
00415 }else{
00416 edm::LogVerbatim ("MainResults") << "else if(fabs(eta)<"<<etabinVals_[etab+1] <<") res = " <<
00417 fResEtaBin[ro][etab]->GetParameter(0) << "+" << fResEtaBin[ro][etab]->GetParameter(1)
00418 << "*exp(-(" <<fResEtaBin[ro][etab]->GetParameter(2) << "*pt));";
00419 }
00420 }else if(nrFilled != 0 && ro == 6){
00421 if(etab == 0){
00422 edm::LogVerbatim ("MainResults") << "if(fabs(eta)<"<<etabinVals_[etab+1] <<") res = " <<
00423 fResEtaBin[ro][etab]->GetParameter(0) << "+" << fResEtaBin[ro][etab]->GetParameter(1)
00424 << "*pt;";
00425 }else{
00426 edm::LogVerbatim ("MainResults") << "else if(fabs(eta)<"<<etabinVals_[etab+1] <<") res = " <<
00427 fResEtaBin[ro][etab]->GetParameter(0) << "+" << fResEtaBin[ro][etab]->GetParameter(1)
00428 << "*pt;";
00429
00430 }
00431 }
00432 }
00433 }
00434 }else if(nrFilled != 0 && objectType_ == "met"){
00435 for(ro=0; ro<8; ro++) {
00436 edm::LogVerbatim ("MainResults") << "-------------------- ";
00437 edm::LogVerbatim ("MainResults") << "\n Resolutions on " << resObsName2[ro] << "\n";
00438 edm::LogVerbatim ("MainResults") << "-------------------- ";
00439 if(nrFilled != 0 && ro != 6) {
00440 edm::LogVerbatim ("MainResults") << "res = " <<
00441 fResEtaBin[ro][0]->GetParameter(0) << "+" << fResEtaBin[ro][0]->GetParameter(1)
00442 << "*exp(-(" <<fResEtaBin[ro][0]->GetParameter(2) << "*pt));";
00443 }else if(nrFilled != 0 && ro == 6){
00444 edm::LogVerbatim ("MainResults") << "res = " <<
00445 fResEtaBin[ro][0]->GetParameter(0) << "+" << fResEtaBin[ro][0]->GetParameter(1) << "*pt;";
00446 }
00447 }
00448 }
00449 }
00450
00451
00452 DEFINE_FWK_MODULE(ResolutionCreator);