00001 #include "RecoParticleFlow/PFClusterTools/interface/PFClusterCalibration.h"
00002 #include "DataFormats/ParticleFlowReco/interface/CalibrationProvenance.h"
00003 #include "DataFormats/ParticleFlowReco/interface/Calibratable.h"
00004
00005 #include <cmath>
00006 #include <cassert>
00007 #include <TBranch.h>
00008 #include <TF1.h>
00009 #include <TTree.h>
00010
00011 using namespace pftools;
00012
00013 void PFClusterCalibration::init() {
00014
00015
00016 correction_ = new TF1("correction",
00017 "((x-[0])/[1])*(x>[4])+((x-[2])/[3])*(x<[4])");
00018 etaCorrection_
00019 = new TF1( "etaCorrection",
00020 "(1-[0]*x-[1]*x*x)*(x<[2])+([3]+[4]*x)*(x>[2]&&x<[5])+(1-[6]*x-[7]*x*x)*(x>[5])");
00021
00022 correction_->FixParameter(0, globalP0_);
00023 correction_->FixParameter(1, globalP1_);
00024 correction_->FixParameter(2, lowEP0_);
00025 correction_->FixParameter(3, lowEP1_);
00026 correction_->FixParameter(4, correctionLowLimit_);
00027
00028
00029
00030
00031
00032 std::string eheb("ecalHcalEcalBarrel");
00033 names_.push_back(eheb);
00034 std::string ehee("ecalHcalEcalEndcap");
00035 names_.push_back(ehee);
00036 std::string ehhb("ecalHcalHcalBarrel");
00037 names_.push_back(ehhb);
00038 std::string ehhe("ecalHcalHcalEndcap");
00039 names_.push_back(ehhe);
00040
00041
00042
00043
00044
00045 const char*
00046 funcString("([0]*[5]*x)*([5]*x<=[1])+([2]+[3]*exp([4]*[5]*x))*([5]*x>[1])");
00047
00048
00049 for (std::vector<std::string>::const_iterator cit = names_.begin(); cit
00050 != names_.end(); ++cit) {
00051 std::string name = *cit;
00052 TF1 func(name.c_str(), funcString);
00053
00054 func.FixParameter(0, 1);
00055 func.FixParameter(1, 0);
00056 func.FixParameter(2, 1);
00057 func.FixParameter(3, 0);
00058 func.FixParameter(4, 0);
00059 func.FixParameter(5, 1);
00060
00061 func.SetMinimum(0);
00062
00063 namesAndFunctions_[name] = func;
00064
00065 }
00066 }
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104 PFClusterCalibration::PFClusterCalibration() :
00105 barrelEndcapEtaDiv_(1.0), ecalOnlyDiv_(0.3), hcalOnlyDiv_(0.5),
00106 doCorrection_(1), allowNegativeEnergy_(0), doEtaCorrection_(1),
00107 maxEToCorrect_(-1.0), correctionLowLimit_(0.), globalP0_(0.0),
00108 globalP1_(1.0), lowEP0_(0.0), lowEP1_(1.0) {
00109
00110 init();
00111
00112
00113
00114
00115 }
00116
00117 PFClusterCalibration::~PFClusterCalibration() {
00118 delete correction_;
00119 delete etaCorrection_;
00120 }
00121
00122 void PFClusterCalibration::setEtaCorrectionParameters(std::vector<double> params) {
00123 if (params.size() != 6) {
00124 std::cout << __PRETTY_FUNCTION__ << ": params is of the wrong length."
00125 << std::endl;
00126 return;
00127 }
00128
00129 unsigned count(0);
00130 for (std::vector<double>::const_iterator dit = params.begin(); dit
00131 != params.end(); ++dit) {
00132
00133 etaCorrection_->FixParameter(count, *dit);
00134 ++count;
00135 }
00136
00137
00138
00139
00140 }
00141
00142 void PFClusterCalibration::setEvolutionParameters(const std::string& sector,
00143 std::vector<double> params) {
00144 TF1* func = &(namesAndFunctions_.find(sector)->second);
00145 unsigned count(0);
00146
00147 for (std::vector<double>::const_iterator dit = params.begin(); dit
00148 != params.end(); ++dit) {
00149 func->FixParameter(count, *dit);
00150
00151 ++count;
00152 }
00153
00154 func->SetMinimum(0);
00155 }
00156
00157 void PFClusterCalibration::setCorrections(const double& lowEP0,
00158 const double& lowEP1, const double& globalP0, const double& globalP1) {
00159
00160 globalP0_ = globalP0;
00161 globalP1_ = globalP1;
00162
00163 lowEP0_ = lowEP0;
00164 lowEP1_ = lowEP1;
00165
00166 correctionLowLimit_ = (lowEP0_ - globalP0_)/(globalP1_ - lowEP1_);
00167
00168 correction_->FixParameter(0, globalP0_);
00169 correction_->FixParameter(1, globalP1_);
00170 correction_->FixParameter(2, lowEP0_);
00171 correction_->FixParameter(3, lowEP1_);
00172 correction_->FixParameter(4, correctionLowLimit_);
00173
00174
00175
00176 }
00177
00178 double PFClusterCalibration::getCalibratedEcalEnergy(const double& ecalE,
00179 const double& hcalE, const double& eta, const double& phi) const {
00180 const TF1* theFunction(0);
00181
00182 if (fabs(eta) < barrelEndcapEtaDiv_) {
00183
00184 theFunction = &(namesAndFunctions_.find("ecalHcalEcalBarrel")->second);
00185 } else {
00186
00187 theFunction = &(namesAndFunctions_.find("ecalHcalEcalEndcap")->second);
00188 }
00189
00190 assert(theFunction != 0);
00191 double totalE(ecalE + hcalE);
00192 double bCoeff = theFunction->Eval(totalE);
00193 return ecalE * bCoeff;
00194 }
00195
00196 double PFClusterCalibration::getCalibratedHcalEnergy(const double& ecalE,
00197 const double& hcalE, const double& eta, const double& phi) const {
00198 const TF1* theFunction(0);
00199
00200 if (fabs(eta) < barrelEndcapEtaDiv_) {
00201
00202 theFunction = &(namesAndFunctions_.find("ecalHcalHcalBarrel")->second);
00203 } else {
00204
00205 theFunction = &(namesAndFunctions_.find("ecalHcalHcalEndcap")->second);
00206 }
00207
00208 double totalE(ecalE + hcalE);
00209 assert(theFunction != 0);
00210 double cCoeff = theFunction->Eval(totalE);
00211 return hcalE * cCoeff;
00212 }
00213
00214 double PFClusterCalibration::getCalibratedEnergy(const double& ecalE,
00215 const double& hcalE, const double& eta, const double& phi) const {
00216 double totalE(ecalE + hcalE);
00217 double answer(totalE);
00218
00219 answer = getCalibratedEcalEnergy(ecalE, hcalE, eta, phi)
00220 + getCalibratedHcalEnergy(ecalE, hcalE, eta, phi);
00221 if (doEtaCorrection_)
00222 answer = answer/etaCorrection_->Eval(eta);
00223
00224 if (maxEToCorrect_> 0 && answer < maxEToCorrect_)
00225 return correction_->Eval(answer);
00226 if (doCorrection_) {
00227 if (maxEToCorrect_> 0 && answer < maxEToCorrect_)
00228 answer = correction_->Eval(answer);
00229 else if (maxEToCorrect_ < 0) {
00230 answer = correction_->Eval(answer);
00231 }
00232 }
00233 if (!allowNegativeEnergy_ && answer < 0)
00234 return 0;
00235 return answer;
00236
00237 }
00238
00239 void PFClusterCalibration::getCalibratedEnergyEmbedAInHcal(double& ecalE,
00240 double& hcalE, const double& eta, const double& phi) const {
00241
00242 double ecalEOld(ecalE);
00243 double hcalEOld(hcalE);
00244
00245 ecalE = getCalibratedEcalEnergy(ecalEOld, hcalEOld, eta, phi);
00246 hcalE = getCalibratedHcalEnergy(ecalEOld, hcalEOld, eta, phi);
00247
00248 double preCorrection(ecalE + hcalE);
00249 if (doEtaCorrection_)
00250 preCorrection = preCorrection/etaCorrection_->Eval(eta);
00251
00252 if (doCorrection_) {
00253 double corrE = correction_->Eval(preCorrection);
00254
00255 double a = corrE - preCorrection;
00256 hcalE += a;
00257 }
00258 if (hcalE < 0 && !allowNegativeEnergy_)
00259 hcalE = 0;
00260
00261 }
00262
00263 void PFClusterCalibration::calibrate(Calibratable& c) {
00264 CalibrationResultWrapper crw;
00265 getCalibrationResultWrapper(c, crw);
00266 c.calibrations_.push_back(crw);
00267
00268 }
00269
00270 void PFClusterCalibration::getCalibrationResultWrapper(const Calibratable& c,
00271 CalibrationResultWrapper& crw) {
00272
00273 crw.ecalEnergy_ = getCalibratedEcalEnergy(c.cluster_energyEcal_,
00274 c.cluster_energyHcal_, fabs(c.cluster_meanEcal_.eta_),
00275 fabs(c.cluster_meanEcal_.phi_));
00276
00277 crw.hcalEnergy_ = getCalibratedHcalEnergy(c.cluster_energyEcal_,
00278 c.cluster_energyHcal_, fabs(c.cluster_meanEcal_.eta_),
00279 fabs(c.cluster_meanEcal_.phi_));
00280
00281 crw.particleEnergy_ = getCalibratedEnergy(c.cluster_energyEcal_,
00282 c.cluster_energyHcal_, fabs(c.cluster_meanEcal_.eta_),
00283 fabs(c.cluster_meanEcal_.phi_));
00284
00285 crw.provenance_ = LINEARCORR;
00286 crw.b_ = crw.ecalEnergy_ / c.cluster_energyEcal_;
00287 crw.c_ = crw.hcalEnergy_ / c.cluster_energyHcal_;
00288 crw.truthEnergy_ = c.sim_energyEvent_;
00289
00290 }
00291
00292 void PFClusterCalibration::calibrateTree(TTree* input) {
00293 std::cout << __PRETTY_FUNCTION__
00294 << ": WARNING! This isn't working properly yet!\n";
00295 TBranch* calibBr = input->GetBranch("Calibratable");
00296 Calibratable* calib_ptr = new Calibratable();
00297 calibBr->SetAddress(&calib_ptr);
00298
00299
00300
00301
00302 std::cout << "Looping over tree's "<< input->GetEntries()
00303 << " entries...\n";
00304 for (unsigned entries(0); entries < 20000; entries++) {
00305 if (entries % 10000== 0)
00306 std::cout << "\tProcessing entry "<< entries << "\n";
00307 input->GetEntry(entries);
00308 calibrate(*calib_ptr);
00309 input->Fill();
00310 }
00311
00312 }
00313
00314 std::ostream& pftools::operator<<(std::ostream& s, const PFClusterCalibration& cc) {
00315 s << "PFClusterCalibration: dump...\n";
00316 s << "barrelEndcapEtaDiv:\t" << cc.barrelEndcapEtaDiv_ << ", ecalOnlyDiv:\t" << cc.ecalOnlyDiv_;
00317 s << ", \nhcalOnlyDiv:\t" << cc.hcalOnlyDiv_ << ", doCorrection:\t" << cc.doCorrection_;
00318 s << ", \nallowNegativeEnergy:\t" << cc.allowNegativeEnergy_;
00319 s << ", \ncorrectionLowLimit:\t" << cc.correctionLowLimit_ << ",parameters:\t" << cc.globalP0_ << ", ";
00320 s << cc.globalP1_ << ", " << cc.lowEP0_ << ", " << cc.lowEP1_;
00321 s << "\ndoEtaCorrection:\t" << cc.doEtaCorrection_;
00322 return s;
00323 }
00324