CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/src/RecoParticleFlow/PFClusterTools/src/CalibCompare.cc

Go to the documentation of this file.
00001 #include "RecoParticleFlow/PFClusterTools/interface/CalibCompare.h"
00002 #include "RecoParticleFlow/PFClusterTools/interface/PFToolsException.h"
00003 #include "RecoParticleFlow/PFClusterTools/interface/TreeUtility.h"
00004 #include "DataFormats/ParticleFlowReco/interface/Calibratable.h"
00005 #include "RecoParticleFlow/PFClusterTools/interface/DetectorElement.h"
00006 #include "RecoParticleFlow/PFClusterTools/interface/Calibrator.h"
00007 #include "RecoParticleFlow/PFClusterTools/interface/LinearCalibrator.h"
00008 #include "RecoParticleFlow/PFClusterTools/interface/SpaceManager.h"
00009 #include "DataFormats/ParticleFlowReco/interface/CalibrationProvenance.h"
00010 #include "RecoParticleFlow/PFClusterTools/interface/Region.h"
00011 #include "RecoParticleFlow/PFClusterTools/interface/IO.h"
00012 
00013 #include <TFile.h>
00014 #include <TTree.h>
00015 #include <TH3F.h>
00016 #include <TF2.h>
00017 #include <TH1F.h>
00018 #include <TGraph.h>
00019 #include <TF1.h>
00020 #include <vector>
00021 #include <TROOT.h>
00022 
00023 using namespace pftools;
00024 
00025 CalibCompare::~CalibCompare() {
00026 }
00027 
00028 CalibCompare::CalibCompare(IO* options) :
00029         withOffset_(false), target_(CLUSTER), options_(options), debug_(0),
00030                         mlpOffset_(0.0), mlpSlope_(1.0) {
00031 
00032         options_->GetOpt("exercises", "withOffset", withOffset_);
00033         options_->GetOpt("exercises", "debug", debug_);
00034 
00035         /* Initialise PFClusterCalibration appropriately. */
00036         double g0, g1, e0, e1;
00037         options_->GetOpt("correction", "globalP0", g0);
00038         options_->GetOpt("correction", "globalP1", g1);
00039         options_->GetOpt("correction", "lowEP0", e0);
00040         options_->GetOpt("correction", "lowEP1", e1);
00041         clusterCalibration_.setCorrections(e0, e1, g0, g1);
00042 
00043         double ecalECut, hcalECut;
00044         options_->GetOpt("evolution", "ecalECut", ecalECut);
00045         options_->GetOpt("evolution", "hcalECut", hcalECut);
00046         clusterCalibration_.setEcalHcalEnergyCuts(ecalECut, hcalECut);
00047 
00048         int allowNegative(0);
00049         options_->GetOpt("correction", "allowNegativeEnergy", allowNegative);
00050         clusterCalibration_.setAllowNegativeEnergy(allowNegative);
00051 
00052         int doCorrection(1);
00053         options_->GetOpt("correction", "doCorrection", doCorrection);
00054         clusterCalibration_.setDoCorrection(doCorrection);
00055 
00056         int doEtaCorrection(0);
00057         options_->GetOpt("correction", "doEtaCorrection", doEtaCorrection);
00058         clusterCalibration_.setDoEtaCorrection(doEtaCorrection);
00059 
00060         double barrelEta;
00061         options_->GetOpt("evolution", "barrelEndcapEtaDiv", barrelEta);
00062         clusterCalibration_.setBarrelBoundary(barrelEta);
00063 
00064         double maxEToCorrect(100.0);
00065         options_->GetOpt("correction", "maxEToCorrect", maxEToCorrect);
00066         clusterCalibration_.setMaxEToCorrect(maxEToCorrect);
00067 
00068         std::vector<std::string>* names = clusterCalibration_.getKnownSectorNames();
00069         for (std::vector<std::string>::iterator i = names->begin(); i
00070                         != names->end(); ++i) {
00071                 std::string sector = *i;
00072                 std::vector<double> params;
00073                 options_->GetOpt("evolution", sector.c_str(), params);
00074                 clusterCalibration_.setEvolutionParameters(sector, params);
00075         }
00076 
00077         std::vector<double> etaParams;
00078         options_->GetOpt("evolution", "etaCorrection", etaParams);
00079         clusterCalibration_.setEtaCorrectionParameters(etaParams);
00080 
00081         std::cout << clusterCalibration_ << "\n";
00082 
00083         options_->GetOpt("correction", "mlpOffset", mlpOffset_);
00084         options_->GetOpt("correction", "mlpSlope", mlpSlope_);
00085 
00086         erlCalibration_.setOffsetAndSlope(mlpOffset_, mlpSlope_);
00087 
00088         if (debug_ > 0)
00089                 std::cout << __PRETTY_FUNCTION__ << ": finished.\n";
00090 
00091 }
00092 
00093 void CalibCompare::calibrateCalibratables(TChain& sourceTree,
00094                 const std::string& exercisefile) {
00095 
00096         if (debug_ > 0) {
00097                 std::cout << "Welcome to " << __PRETTY_FUNCTION__ << "\n";
00098                 std::cout << "Opening TTree...\n";
00099         }
00100 
00101         TreeUtility tu;
00102         std::vector<Calibratable> calibVec;
00103 
00104         tu.getCalibratablesFromRootFile(sourceTree, calibVec);
00105 
00106         std::cout << "Moving on... " << std::endl;
00107         TFile* exercises = new TFile(exercisefile.c_str(), "recreate");
00108         TTree tree("CalibratedParticles", "");
00109         Calibratable* calibrated = new Calibratable();
00110         tree.Branch("Calibratable", "pftools::Calibratable", &calibrated, 32000, 2);
00111 
00112         evaluateCalibrations(tree, calibrated, calibVec);
00113 
00114         //save results
00115         std::cout << "Writing output tree...\n";
00116         tree.Write();
00117         //gaussianFits(*exercises, calibVec);
00118         exercises->Write();
00119         exercises->Close();
00120         std::cout << "Done." << std::endl;
00121         delete calibrated;
00122 
00123 }
00124 
00125 void CalibCompare::evaluateCalibrations(TTree& tree, Calibratable* calibrated,
00126                 const std::vector<Calibratable>& calibVec) {
00127 
00128         unsigned count(0);
00129         for (std::vector<Calibratable>::const_iterator zit = calibVec.begin(); zit
00130                         != calibVec.end(); ++zit) {
00131 
00132                 const Calibratable& calib = *zit;
00133 
00134                 calibrated->reset();
00135 
00136                 CalibrationResultWrapper crwPre;
00137                 crwPre.ecalEnergy_ = calib.cluster_energyEcal_;
00138                 crwPre.hcalEnergy_ = calib.cluster_energyHcal_;
00139                 crwPre.particleEnergy_ = calib.cluster_energyEcal_
00140                                 + calib.cluster_energyHcal_;
00141                 crwPre.truthEnergy_ = calib.sim_energyEvent_;
00142                 crwPre.provenance_ = UNCALIBRATED;
00143                 crwPre.targetFuncContrib_ = 0;
00144                 crwPre.target_ = target_;
00145                 crwPre.compute();
00146                 calibrated->calibrations_.push_back(crwPre);
00147 
00148                 CalibrationResultWrapper crwErl;
00149 
00150                 crwErl.particleEnergy_ = erlCalibration_.evaluate(crwPre.ecalEnergy_,
00151                                 crwPre.hcalEnergy_, calib.cluster_numEcal_,
00152                                 calib.cluster_numHcal_, fabs(calib.cluster_meanEcal_.eta_)
00153                                                 / 2.0, crwPre.ecalEnergy_ / (crwPre.particleEnergy_),
00154                                 (calib.cluster_meanEcal_.phi_ + 3.14) / 6.3);
00155                 crwErl.ecalEnergy_ = crwErl.particleEnergy_
00156                                 * erlCalibration_.ecalFraction(crwPre.ecalEnergy_,
00157                                                 crwPre.hcalEnergy_, calib.cluster_numEcal_,
00158                                                 calib.cluster_numHcal_, fabs(
00159                                                                 calib.cluster_meanEcal_.eta_) / 2.0,
00160                                                 crwPre.ecalEnergy_ / (crwPre.particleEnergy_),
00161                                                 (calib.cluster_meanEcal_.phi_ + 3.14) / 6.3);
00162 
00163                 crwErl.hcalEnergy_ = crwErl.particleEnergy_ - crwErl.ecalEnergy_;
00164                 crwErl.b_ = crwErl.ecalEnergy_ / crwPre.ecalEnergy_;
00165                 crwErl.c_ = crwErl.hcalEnergy_ / crwPre.hcalEnergy_;
00166 
00167                 crwErl.truthEnergy_ = calib.sim_energyEvent_;
00168                 crwErl.provenance_ = BAYESIAN;
00169                 crwErl.target_ = target_;
00170                 crwErl.compute();
00171                 calibrated->calibrations_.push_back(crwErl);
00172 
00173                 CalibrationResultWrapper crwCorr;
00174 
00175                 crwCorr.ecalEnergy_ = clusterCalibration_.getCalibratedEcalEnergy(
00176                                 crwPre.ecalEnergy_, crwPre.hcalEnergy_,
00177                                 calib.cluster_meanEcal_.eta_, calib.cluster_meanEcal_.phi_);
00178                 crwCorr.hcalEnergy_ = clusterCalibration_.getCalibratedHcalEnergy(
00179                                 crwPre.ecalEnergy_, crwPre.hcalEnergy_,
00180                                 calib.cluster_meanHcal_.eta_, calib.cluster_meanHcal_.phi_);
00181                 crwCorr.particleEnergy_ = clusterCalibration_.getCalibratedEnergy(
00182                                 crwPre.ecalEnergy_, crwPre.hcalEnergy_, calib.sim_etaEcal_,
00183                                 calib.sim_phiEcal_);
00184 
00185                 crwCorr.truthEnergy_ = calib.sim_energyEvent_;
00186                 crwCorr.provenance_ = LINEAR;
00187                 crwCorr.targetFuncContrib_ = 0;;
00188                 crwCorr.target_ = target_;
00189                 crwCorr.compute();
00190                 calibrated->calibrations_.push_back(crwCorr);
00191 
00192                 calibrated->recompute();
00193 
00194                 tree.Fill();
00195 
00196                 ++count;
00197         }
00198 
00199 }
00200