CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_4/src/RecoParticleFlow/PFClusterTools/src/Exercises3.cc

Go to the documentation of this file.
00001 #include "RecoParticleFlow/PFClusterTools/interface/Exercises3.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 void resetElement3(DetectorElementPtr de) {
00026         de->setCalib(1.0);
00027 }
00028 
00029 Exercises3::~Exercises3() {
00030         calibResultsFile_.close();
00031 }
00032 
00033 Exercises3::Exercises3(IO* options) :
00034         withOffset_(false), target_(CLUSTER), threshold_(30), options_(options),
00035                         debug_(0) {
00036 
00037         options_->GetOpt("exercises", "withOffset", withOffset_);
00038         options_->GetOpt("exercises", "threshold", threshold_);
00039         options_->GetOpt("exercises", "debug", debug_);
00040 
00041         /* Initialise PFClusterCalibration appropriately. */
00042         double g0, g1, e0, e1;
00043         options_->GetOpt("correction", "globalP0", g0);
00044         options_->GetOpt("correction", "globalP1", g1);
00045         options_->GetOpt("correction", "lowEP0", e0);
00046         options_->GetOpt("correction", "lowEP1", e1);
00047         clusterCalibration_.setCorrections(e0, e1, g0, g1);
00048 
00049         double ecalECut, hcalECut;
00050         options_->GetOpt("evolution", "ecalECut", ecalECut);
00051         options_->GetOpt("evolution", "hcalECut", hcalECut);
00052         clusterCalibration_.setEcalHcalEnergyCuts(ecalECut, hcalECut);
00053 
00054         int allowNegative(0);
00055         options_->GetOpt("correction", "allowNegativeEnergy", allowNegative);
00056         clusterCalibration_.setAllowNegativeEnergy(allowNegative);
00057 
00058         int doCorrection(1);
00059         options_->GetOpt("correction", "doCorrection", doCorrection);
00060         clusterCalibration_.setDoCorrection(doCorrection);
00061 
00062         int doEtaCorrection(0);
00063         options_->GetOpt("correction", "doEtaCorrection", doEtaCorrection);
00064         clusterCalibration_.setDoEtaCorrection(doEtaCorrection);
00065 
00066         double barrelEta;
00067         options_->GetOpt("evolution", "barrelEndcapEtaDiv", barrelEta);
00068         clusterCalibration_.setBarrelBoundary(barrelEta);
00069 
00070         double maxEToCorrect(100.0);
00071         options_->GetOpt("correction", "maxEToCorrect", maxEToCorrect);
00072         clusterCalibration_.setMaxEToCorrect(maxEToCorrect);
00073 
00074         std::vector<std::string>* names = clusterCalibration_.getKnownSectorNames();
00075         for (std::vector<std::string>::iterator i = names->begin(); i
00076                         != names->end(); ++i) {
00077                 std::string sector = *i;
00078                 std::vector<double> params;
00079                 options_->GetOpt("evolution", sector.c_str(), params);
00080                 clusterCalibration_.setEvolutionParameters(sector, params);
00081         }
00082 
00083         std::vector<double> etaParams;
00084         options_->GetOpt("evolution", "etaCorrection", etaParams);
00085         clusterCalibration_.setEtaCorrectionParameters(etaParams);
00086 
00087         std::cout << clusterCalibration_ << "\n";
00088 
00089         std::string outputFileName;
00090         options_->GetOpt("results", "calibParamOutput", outputFileName);
00091         calibResultsFile_.open(outputFileName.c_str());
00092         calibResultsFile_ << "//Hello from your friendly PFClusterTools!\n";
00093         if (debug_ > 0)
00094                 std::cout << __PRETTY_FUNCTION__ << ": finished.\n";
00095 
00096 }
00097 
00098 void Exercises3::calibrateCalibratables(TChain& sourceTree,
00099                 const std::string& exercisefile) {
00100 
00101         if (debug_ > 0) {
00102                 std::cout << "Welcome to "<< __PRETTY_FUNCTION__ << "\n";
00103                 std::cout << "Opening TTree...\n";
00104         }
00105 //      //open tfile
00106 //      TFile* source = new TFile(sourcefile.c_str());
00107 //      if (source == 0) {
00108 //              std::string desc("Couldn't open file ");
00109 //              desc += sourcefile;
00110 //              PFToolsException e(desc);
00111 //              throw e;
00112 //      }
00113 //      if (debug_ > 0)
00114 //              std::cout << "Extracting calibratables...\n";
00115         //use tree utility to extract calibratables
00116         TreeUtility tu;
00117         std::vector<Calibratable> calibVec;
00118 
00119         tu.getCalibratablesFromRootFile(sourceTree, calibVec);
00120         if (debug_ > 0)
00121                 std::cout << "Got a vector of calibratables of size "<< calibVec.size()
00122                                 << "\n";
00123         //initialise detector elements
00124         DetectorElementPtr ecal(new DetectorElement(ECAL, 1.0));
00125         DetectorElementPtr hcal(new DetectorElement(HCAL, 1.0));
00126         DetectorElementPtr offset(new DetectorElement(OFFSET, 1.0));
00127 
00128         //convert calibratables to particle deposits
00129         std::vector<ParticleDepositPtr> pdVec;
00130         tu.convertCalibratablesToParticleDeposits(calibVec, pdVec, target_, offset,
00131                         ecal, hcal, withOffset_);
00132         //source->Close();
00133 
00134         TFile* exercises = new TFile(exercisefile.c_str(), "recreate");
00135         TH1F droppedParticles("droppedParticles", "droppedParticles", 100000, 0,
00136                         100000);
00137         if (debug_ > 0)
00138                 std::cout << "Particle deposit vec has "<< pdVec.size() << " entries\n";
00139 
00140         //calibrate
00141         if (debug_ > 1)
00142                 std::cout << "Creating calibrator clones and space managers\n";
00143         boost::shared_ptr<Calibrator> linCal(new LinearCalibrator());
00144 
00145         //Tell the calibrator which detector elements should be calibrated
00146         if (withOffset_) {
00147                 linCal->addDetectorElement(offset);
00148 
00149         }
00150         linCal->addDetectorElement(ecal);
00151         linCal->addDetectorElement(hcal);
00152 
00153         double barrelEta;
00154         options_->GetOpt("evolution", "barrelEndcapEtaDiv", barrelEta);
00155         boost::shared_ptr<SpaceManager> sm(new SpaceManager("ecalAndHcal"));
00156         sm->setBarrelLimit(barrelEta);
00157         sm->createCalibrators(*linCal);
00158 
00159         if (debug_ > 1)
00160                 std::cout << "Initialised SpaceManager and calibrators.\n";
00161         elements_.clear();
00162         if (withOffset_)
00163                 elements_.push_back(offset);
00164         elements_.push_back(ecal);
00165         elements_.push_back(hcal);
00166 
00167         //Initialise calibrators with particles
00168         int dropped(0);
00169 
00170         double eCut(0.0);
00171         double hCut(0.0);
00172         options_->GetOpt("evolution", "ecalECut", eCut);
00173         options_->GetOpt("evolution", "hcalECut", hCut);
00174         if (debug_ > 0)
00175                 std::cout << "Using a ECAL and HCAL energy cuts of "<< eCut << " and "
00176                                 << hCut << " GeV\n";
00177         if (debug_ > 1)
00178                 std::cout << "Assigning particles to space managers and calibrators.\n";
00179 
00180         //This is just a convenience plot to check on the hcal
00181         for (std::vector<ParticleDepositPtr>::const_iterator cit = pdVec.begin(); cit
00182                         != pdVec.end(); ++cit) {
00183                 ParticleDepositPtr pd = *cit;
00184                 CalibratorPtr c = sm->findCalibrator(pd->getEta(), pd->getPhi(),
00185                                 pd->getTruthEnergy());
00186                 //std::cout << *pd << "\n";
00187                 if (c == 0) {
00188                         std::cout << "Couldn't find calibrator for particle?!\n";
00189                         std::cout << "\t"<< *pd << "\n";
00190 
00191                         dropped++;
00192                 } else {
00193                         c->addParticleDeposit(pd);
00194                 }
00195 
00196         }
00197 
00198         if (debug_ > 1)
00199                 std::cout << "Dropped "<< dropped << " particles.\n";
00200 
00201         /* Done assignments, now calibrate */
00202         if (debug_ > 1)
00203                 std::cout
00204                                 << "Assignments complete, starting calibration and analysis.\n";
00205 
00206         //calibrate
00207         std::map<SpaceVoxelPtr, CalibratorPtr> * smCalibrators =
00208                         sm->getCalibrators();
00209 
00210         TTree tree("CalibratedParticles", "");
00211         Calibratable* calibrated = new Calibratable();
00212         tree.Branch("Calibratable", "pftools::Calibratable", &calibrated, 32000, 2);
00213         if (debug_ > 1)
00214                 std::cout << "Initialised tree.\n";
00215 
00216         /* ECAL and HCAL */
00217         std::cout << "*** Performance for ECAL + HCAL calibration ***\n";
00218         getCalibrations(sm);
00219         exercises->cd("/");
00220         exercises->mkdir("ecalAndHcal");
00221         exercises->cd("/ecalAndHcal");
00222         evaluateSpaceManager(sm, elements_);
00223         for (std::map<SpaceVoxelPtr, CalibratorPtr>::iterator it =
00224                         smCalibrators->begin(); it != smCalibrators->end(); ++it) {
00225                 SpaceVoxelPtr sv = (*it).first;
00226                 CalibratorPtr c = (*it).second;
00227                 std::for_each(elements_.begin(), elements_.end(), resetElement3);
00228                 evaluateCalibrator(sm, c, tree, calibrated, ecal, hcal, offset, LINEAR,
00229                                 LINEARCORR);
00230 
00231                 std::for_each(elements_.begin(), elements_.end(), resetElement3);
00232         }
00233         sm->printCalibrations(std::cout);
00234 
00235         exercises->cd("/");
00236 
00237         //save results
00238         std::cout << "Writing output tree...\n";
00239         tree.Write();
00240         droppedParticles.Write();
00241         //gaussianFits(*exercises, calibVec);
00242         exercises->Write();
00243         exercises->Close();
00244         std::cout << "Done."<< std::endl;
00245 
00246 }
00247 
00248 void Exercises3::getCalibrations(SpaceManagerPtr s) {
00249 
00250         std::map<SpaceVoxelPtr, CalibratorPtr>* smCalibrators = s->getCalibrators();
00251 
00252         for (std::map<SpaceVoxelPtr, CalibratorPtr>::iterator it =
00253                         smCalibrators->begin(); it != smCalibrators->end(); ++it) {
00254                 CalibratorPtr c= (*it).second;
00255                 std::for_each(elements_.begin(), elements_.end(), resetElement3);
00256                 if (c->hasParticles() > static_cast<int>(threshold_)) {
00257                         std::map<DetectorElementPtr, double> calibs =
00258                                         c->getCalibrationCoefficients();
00259                         s->assignCalibration(c, calibs);
00260                 }
00261         }
00262 }
00263 
00264 void Exercises3::evaluateSpaceManager(SpaceManagerPtr s,
00265                 std::vector<DetectorElementPtr> detEls) {
00266 
00267         int autoFit(0);
00268         options_->GetOpt("evolution", "autoFit", autoFit);
00269         std::cout << "AutoFit option = "<< autoFit << "\n";
00270 
00271         std::vector<double> ecalBarrel;
00272         std::vector<double> ecalEndcap;
00273         std::vector<double> hcalBarrel;
00274         std::vector<double> hcalEndcap;
00275 
00276         double minE, maxE;
00277         options_->GetOpt("evolution", "evolutionFunctionMinE", minE);
00278         options_->GetOpt("evolution", "evolutionFunctionMaxE", maxE);
00279 
00280         int basePlots(0);
00281         options_->GetOpt("evolution", "basePlots", basePlots);
00282 
00283         int useTruth(1);
00284         options_->GetOpt("evolution", "basePlotsUseTruth", useTruth);
00285 
00286         if (debug_ > 1&& basePlots > 0)
00287                 std::cout << "Option to generate evolution plots invoked.\n";
00288 
00289         if (autoFit == 0) {
00290                 std::cout << "Fixing parameters for evolution functions.\n";
00291 
00292                 options_->GetOpt("evolution", "ecalHcalEcalBarrel", ecalBarrel);
00293                 options_->GetOpt("evolution", "ecalHcalEcalEndcap", ecalEndcap);
00294                 options_->GetOpt("evolution", "ecalHcalHcalBarrel", hcalBarrel);
00295                 options_->GetOpt("evolution", "ecalHcalHcalEndcap", hcalEndcap);
00296                 assert(ecalBarrel.size() == 6 && ecalEndcap.size() == 6);
00297                 assert(hcalBarrel.size() == 6 && hcalEndcap.size() == 6);
00298 
00299                 for (std::vector<DetectorElementPtr>::iterator i = detEls.begin(); i
00300                                 != detEls.end(); ++i) {
00301                         DetectorElementPtr d = *i;
00302                         std::cout << "Fixing evolution for "<< *d << "\n";
00303                         int mode(0);
00304                         options_->GetOpt("spaceManager", "interpolationMode", mode);
00305 
00306                         std::string name("Func_");
00307                         name.append(DetElNames[d->getType()]);
00308                         name.append("_");
00309 
00310                         /* Fitting for barrel */
00311                         std::string barrelName(name);
00312                         barrelName.append(RegionNames[BARREL_POS]);
00313                         std::cout << "\tFixing "<< RegionNames[BARREL_POS]<< "\n";
00314                         TF1
00315                                         fBarrel(barrelName.c_str(),
00316                                                         "([0]*[5]*x)*([5]*x<[1])+([2]+[3]*exp([4]*[5]*x))*([5]*x>[1])");
00317 
00318                         if (d->getType() == ECAL) {
00319                                 unsigned count(0);
00320                                 for (std::vector<double>::const_iterator dit =
00321                                                 ecalBarrel.begin(); dit!= ecalBarrel.end(); ++dit) {
00322                                         fBarrel.FixParameter(count, *dit);
00323                                         ++count;
00324                                 }
00325 
00326                         }
00327                         if (d->getType() == HCAL) {
00328                                 unsigned count(0);
00329                                 for (std::vector<double>::const_iterator dit =
00330                                                 hcalBarrel.begin(); dit!= hcalBarrel.end(); ++dit) {
00331                                         fBarrel.FixParameter(count, *dit);
00332                                         ++count;
00333                                 }
00334 
00335                         }
00336                         if (useTruth)
00337                                 fBarrel.FixParameter(5, 1.0);
00338 
00339                         fBarrel.SetMinimum(0);
00340                         s->addEvolution(d, BARREL_POS, fBarrel);
00341 
00342                         if (basePlots > 0) {
00343                                 TH1* slices = s->extractEvolution(d, BARREL_POS, fBarrel,
00344                                                 useTruth);
00345                                 slices->Write();
00346                         }
00347                         fBarrel.Write();
00348 
00349                         /* Fitting for endcap */
00350                         std::string endcapName(name);
00351                         endcapName.append(RegionNames[ENDCAP_POS]);
00352                         std::cout << "\tFixing "<< RegionNames[ENDCAP_POS]<< "\n";
00353                         TF1
00354                                         fEndcap(endcapName.c_str(),
00355                                                         "([0]*[5]*x)*([5]*x<[1])+([2]+[3]*exp([4]*[5]*x))*([5]*x>[1])");
00356                         
00357 
00358                         if (d->getType() == ECAL) {
00359                                 unsigned count(0);
00360                                 for (std::vector<double>::const_iterator dit =
00361                                                 ecalEndcap.begin(); dit!= ecalEndcap.end(); ++dit) {
00362                                         fEndcap.FixParameter(count, *dit);
00363                                         ++count;
00364                                 }
00365 
00366                         }
00367                         if (d->getType() == HCAL) {
00368                                 unsigned count(0);
00369                                 for (std::vector<double>::const_iterator dit =
00370                                                 hcalEndcap.begin(); dit!= hcalEndcap.end(); ++dit) {
00371                                         fEndcap.FixParameter(count, *dit);
00372                                         ++count;
00373                                 }
00374 
00375                         }
00376                         if (useTruth)
00377                                 fEndcap.FixParameter(5, 1.0);
00378 
00379                         fEndcap.SetMinimum(0);
00380                         s->addEvolution(d, ENDCAP_POS, fEndcap);
00381                         if (basePlots > 0) {
00382                                 TH1* slices = s->extractEvolution(d, ENDCAP_POS, fEndcap,
00383                                                 useTruth);
00384                                 slices->Write();
00385                         }
00386                         fEndcap.Write();
00387                 }
00388 
00389         }
00390 
00391 }
00392 
00393 void Exercises3::evaluateCalibrator(SpaceManagerPtr s, CalibratorPtr c,
00394                 TTree& tree, Calibratable* calibrated, DetectorElementPtr ecal,
00395                 DetectorElementPtr hcal, DetectorElementPtr offset,
00396                 CalibrationProvenance cp, CalibrationProvenance cpCorr) {
00397 
00398         if (c->hasParticles() > static_cast<int>(threshold_)) {
00399                 std::map<DetectorElementPtr, double> calibs = s->getCalibration(c);
00400 
00401                 std::vector<ParticleDepositPtr> csParticles = c->getParticles();
00402                 unsigned count(0);
00403                 int mode(0);
00404                 options_->GetOpt("spaceManager", "interpolationMode", mode);
00405                 if (debug_ > 1) {
00406                         std::cout << "Using interpolation mode " << mode << "\n";
00407                 }
00408 
00409                 for (std::vector<ParticleDepositPtr>::iterator zit = csParticles.begin(); zit
00410                                 != csParticles.end(); ++zit) {
00411                         ParticleDepositPtr pd = *zit;
00412                         calibrated->reset();
00413                         calibrated->rechits_meanEcal_.energy_ = pd->getRecEnergy(ecal);
00414                         calibrated->rechits_meanHcal_.energy_ = pd->getRecEnergy(hcal);
00415                         calibrated->sim_energyEvent_ = pd->getTruthEnergy();
00416                         calibrated->sim_etaEcal_ = pd->getEta();
00417 
00418                         for (std::map<DetectorElementPtr, double>::iterator deit =
00419                                         calibs.begin(); deit != calibs.end(); ++deit) {
00420                                 DetectorElementPtr de = (*deit).first;
00421                                 de->setCalib(1.0);
00422                         }
00423 
00424                         CalibrationResultWrapper crwPre;
00425                         crwPre.ecalEnergy_ = pd->getRecEnergy(ecal);
00426                         crwPre.hcalEnergy_ = pd->getRecEnergy(hcal);
00427                         crwPre.particleEnergy_ = pd->getRecEnergy();
00428                         crwPre.truthEnergy_ = pd->getTruthEnergy();
00429                         crwPre.provenance_ = UNCALIBRATED;
00430                         crwPre.targetFuncContrib_ = pd->getTargetFunctionContrib();
00431                         crwPre.target_ = target_;
00432                         crwPre.compute();
00433                         calibrated->calibrations_.push_back(crwPre);
00434 
00435                         double tempEnergy = pd->getRecEnergy();
00436                         //evaluate calibration
00437                         for (std::map<DetectorElementPtr, double>::iterator deit =
00438                                         calibs.begin(); deit != calibs.end(); ++deit) {
00439                                 DetectorElementPtr de = (*deit).first;
00440 
00441                                 if (mode == 1)
00442                                         de->setCalib(s->interpolateCoefficient(de,
00443                                                         pd->getTruthEnergy(), pd->getEta(), pd->getPhi()));
00444                                 else if (mode == 2|| mode == 3|| mode == 4)
00445                                         de->setCalib(s->evolveCoefficient(de, tempEnergy,
00446                                                         pd->getEta(), pd->getPhi()));
00447                                 else
00448                                         de->setCalib((*deit).second);
00449                         }
00450                         if (debug_ > 8) {
00451                                 std::cout << "POST ECAL HCAL coeff: " << ecal->getCalib() << ", " << hcal->getCalib() << "\n";
00452                         }
00453 
00454                         CalibrationResultWrapper crwPos;
00455                         crwPos.ecalEnergy_ = pd->getRecEnergy(ecal);
00456                         crwPos.hcalEnergy_ = pd->getRecEnergy(hcal);
00457                         crwPos.b_ = ecal->getCalib();
00458                         crwPos.c_ = hcal->getCalib();
00459                         crwPos.particleEnergy_ = pd->getRecEnergy();
00460                         crwPos.truthEnergy_ = pd->getTruthEnergy();
00461                         crwPos.provenance_ = cp;
00462                         crwPos.compute();
00463 
00464                         crwPos.targetFuncContrib_ = pd->getTargetFunctionContrib();
00465                         crwPos.target_ = target_;
00466                         calibrated->calibrations_.push_back(crwPos);
00467 
00468                         //same again, but applying correction
00469                         if (cpCorr != NONE) {
00470                                 CalibrationResultWrapper crwCorr;
00471 
00472                                 crwCorr.ecalEnergy_
00473                                                 = clusterCalibration_.getCalibratedEcalEnergy(
00474                                                                 crwPre.ecalEnergy_, crwPre.hcalEnergy_,
00475                                                                 pd->getEta(), pd->getPhi());
00476                                 crwCorr.hcalEnergy_
00477                                                 = clusterCalibration_.getCalibratedHcalEnergy(
00478                                                                 crwPre.ecalEnergy_, crwPre.hcalEnergy_,
00479                                                                 pd->getEta(), pd->getPhi());
00480                                 if (debug_ > 8) {
00481                                         if(crwPre.ecalEnergy_  > 0 && crwPre.hcalEnergy_ >0)
00482                                         std::cout << "CORR ECAL HCAL coeff: " << crwCorr.ecalEnergy_ / crwPre.ecalEnergy_  << ", " << crwCorr.hcalEnergy_/ crwPre.hcalEnergy_ << "\n\n";
00483                                 }
00484 
00485                                 crwCorr.particleEnergy_
00486                                                 = clusterCalibration_.getCalibratedEnergy(
00487                                                                 crwPre.ecalEnergy_, crwPre.hcalEnergy_,
00488                                                                 pd->getEta(), pd->getPhi());
00489 
00490                                 crwCorr.b_ = ecal->getCalib();
00491                                 crwCorr.c_ = hcal->getCalib();
00492 
00493                                 crwCorr.truthEnergy_ = pd->getTruthEnergy();
00494                                 crwCorr.provenance_ = cpCorr;
00495                                 crwCorr.targetFuncContrib_ = pd->getTargetFunctionContrib();
00496                                 crwCorr.target_ = target_;
00497                                 crwCorr.compute();
00498                                 calibrated->calibrations_.push_back(crwCorr);
00499 
00500                                 crwPos.targetFuncContrib_ = pd->getTargetFunctionContrib();
00501                                 crwPos.target_ = target_;
00502                                 calibrated->calibrations_.push_back(crwPos);
00503                         }
00504                         tree.Fill();
00505                         ++count;
00506 
00507                 }
00508         }
00509 }
00510