CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch2/src/RecoParticleFlow/PFClusterTools/plugins/CalibratableTest.cc

Go to the documentation of this file.
00001 #include "RecoParticleFlow/PFClusterTools/interface/CalibratableTest.h"
00002 #include "DataFormats/ParticleFlowReco/interface/PFSimParticle.h"
00003 #include "DataFormats/ParticleFlowReco/interface/PFSimParticleFwd.h"
00004 #include "DataFormats/ParticleFlowReco/interface/PFRecHit.h"
00005 #include "DataFormats/ParticleFlowReco/interface/PFRecHitFwd.h"
00006 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
00007 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
00008 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementFwd.h"
00009 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementCluster.h"
00010 #include "DataFormats/ParticleFlowReco/interface/PFBlock.h"
00011 #include "DataFormats/ParticleFlowReco/interface/PFRecTrackFwd.h"
00012 #include "DataFormats/ParticleFlowReco/interface/PFTrack.h"
00013 #include "DataFormats/ParticleFlowReco/interface/PFNuclearInteraction.h"
00014 #include "DataFormats/ParticleFlowReco/interface/PFClusterFwd.h"
00015 #include "DataFormats/ParticleFlowReco/interface/PFCluster.h"
00016 #include "DataFormats/ParticleFlowReco/interface/PFRecHitFraction.h"
00017 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00018 #include "DataFormats/Math/interface/LorentzVector.h"
00019 
00020 #include "FWCore/PluginManager/interface/ModuleDef.h"
00021 #include "FWCore/Framework/interface/InputSourceMacros.h"
00022 #include "FWCore/Framework/interface/MakerMacros.h"
00023 
00024 using namespace edm;
00025 using namespace pftools;
00026 using namespace reco;
00027 
00028 CalibratableTest::CalibratableTest(const edm::ParameterSet& parameters) :
00029         debug_(0), nParticleWrites_(0), nParticleFails_(0), nEventWrites_(0), nEventFails_(0),
00030                          deltaRCandToSim_(0.4) {
00031 
00032         std::cout << __PRETTY_FUNCTION__ << std::endl;
00033         
00034         /* This procedure is GENERIC to storing any dictionary enable class in a ROOT tree. */
00035         tree_ = fileservice_->make<TTree>("CalibratableTest", "");
00036         calib_ = new Calibratable();
00037         tree_->Branch("Calibratable", "pftools::Calibratable", &calib_, 32000, 2);
00038         
00039         inputTagCandidates_= parameters.getParameter<InputTag>("PFCandidates");
00040         inputTagSimParticles_= parameters.getParameter<InputTag>("PFSimParticles");
00041         inputTagClustersEcal_= parameters.getParameter<InputTag>("PFClustersEcal");
00042         inputTagClustersHcal_= parameters.getParameter<InputTag>("PFClustersHcal");
00043         deltaRCandToSim_ = parameters.getParameter<double>("deltaRCandToSim");
00044         debug_= parameters.getParameter<int>("debug");
00045 }
00046 
00047 CalibratableTest::~CalibratableTest() {
00048 
00049 }
00050 
00051 void CalibratableTest::beginJob() {
00052         if (debug_ > 1)
00053                 std::cout << __PRETTY_FUNCTION__ << "\n";
00054 
00055 }
00056 
00057 void CalibratableTest::analyze(const edm::Event& event,
00058                 const edm::EventSetup& iSetup) {
00059         if (debug_ > 1)
00060                 std::cout << __PRETTY_FUNCTION__ << "\n";
00061 
00062         //Extract new collection references
00063         pfCandidates_ = new Handle<PFCandidateCollection>;
00064         simParticles_ = new Handle<PFSimParticleCollection>;
00065         clustersEcal_ = new Handle<PFClusterCollection>;
00066         clustersHcal_ = new Handle<PFClusterCollection>;
00067 
00068         getCollection(*pfCandidates_, inputTagCandidates_, event);
00069         getCollection(*simParticles_, inputTagSimParticles_, event);
00070         getCollection(*clustersEcal_, inputTagClustersEcal_, event);
00071         getCollection(*clustersHcal_, inputTagClustersHcal_, event);
00072 
00073         //Reset calibratable branch
00074         thisEventPasses_ = true;
00075         thisParticlePasses_ = true;
00076         calib_->reset();
00077 
00078         if (debug_ > 1)
00079                 std::cout << "\tStarting event..."<< std::endl;
00080 
00081         PFSimParticleCollection sims = **simParticles_;
00082         PFCandidateCollection candidates = **pfCandidates_;
00083         PFClusterCollection clustersEcal = **clustersEcal_;
00084         PFClusterCollection clustersHcal = **clustersHcal_;
00085 
00086         if (sims.size() == 0) {
00087                 std::cout << "\tNo sim particles found!" << std::endl;
00088                 thisEventPasses_ = false;
00089         }
00090 
00091         //Find primary pions in the event
00092         std::vector<unsigned> primarySims = findPrimarySimParticles(sims);
00093         if (debug_) {
00094                 std::cout << "\tFound "<< primarySims.size()
00095                                 << " primary sim particles, "<< (**pfCandidates_).size() << " pfCandidates\n";
00096         }
00097         for (std::vector<unsigned>::const_iterator cit = primarySims.begin(); cit
00098                         != primarySims.end(); ++cit) {
00099                 //There will be one write to the tree for each pion found.
00100                 if (debug_ > 1)
00101                         std::cout << "\t**Starting particle...**\n";
00102                 const PFSimParticle& sim = sims[*cit];
00103                 //One sim per calib =>
00104                 calib_->sim_numEvent_ = 1;
00105                 calib_->sim_isMC_ = true;
00106                 calib_->sim_energyEvent_ = sim.trajectoryPoint(PFTrajectoryPoint::ClosestApproach).momentum().E();
00107                 calib_->sim_phi_ = sim.trajectoryPoint(PFTrajectoryPoint::ClosestApproach).momentum().Phi();
00108                 calib_->sim_eta_ = sim.trajectoryPoint(PFTrajectoryPoint::ClosestApproach).momentum().Eta();
00109 
00110                 if (sim.nTrajectoryPoints() > PFTrajectoryPoint::ECALEntrance) {
00111                         calib_->sim_etaEcal_ = sim.trajectoryPoint(PFTrajectoryPoint::ECALEntrance).positionREP().Eta();
00112                         calib_->sim_phiEcal_ = sim.trajectoryPoint(PFTrajectoryPoint::ECALEntrance).positionREP().Phi();
00113                 }
00114                 if (sim.nTrajectoryPoints() > PFTrajectoryPoint::HCALEntrance) {
00115                         calib_->sim_etaHcal_ = sim.trajectoryPoint(PFTrajectoryPoint::HCALEntrance).positionREP().Eta();
00116                         calib_->sim_phiHcal_ = sim.trajectoryPoint(PFTrajectoryPoint::HCALEntrance).positionREP().Phi();
00117                 }
00118 
00119                 // Find candidates near this sim particle
00120                 std::vector<unsigned> matchingCands = findCandidatesInDeltaR(sim,
00121                                 candidates, deltaRCandToSim_);
00122                 if (debug_ > 3)
00123                         std::cout << "\t\tFound candidates near sim, found "
00124                                         << matchingCands.size()<< " of them.\n";
00125                 if (matchingCands.size() == 0)
00126                         thisParticlePasses_ = false;
00127                 for (std::vector<unsigned>::const_iterator mcIt = matchingCands.begin(); mcIt
00128                                 != matchingCands.end(); ++mcIt) {
00129                         const PFCandidate& theCand = candidates[*mcIt];
00130                         extractCandidate(theCand);
00131                 }
00132                 //Finally,
00133                 fillTreeAndReset();
00134                 
00135         }
00136 
00137         delete pfCandidates_;
00138         delete simParticles_;
00139         delete clustersEcal_;
00140         delete clustersHcal_;
00141 
00142         if (thisEventPasses_)
00143                 ++nEventWrites_;
00144         else
00145                 ++nEventFails_;
00146 
00147 }
00148 
00149 std::vector<unsigned> CalibratableTest::findPrimarySimParticles(
00150                 const std::vector<PFSimParticle>& sims) {
00151         std::vector<unsigned> answers;
00152         unsigned index(0);
00153         for (std::vector<PFSimParticle>::const_iterator cit = sims.begin(); cit
00154                         != sims.end(); ++cit) {
00155                 PFSimParticle theSim = *cit;
00156                 //TODO: what about rejected events?
00157                 if (theSim.motherId() >= 0)
00158                         continue;
00159                 int particleId = abs(theSim.pdgCode());
00160                 if (particleId != 211)
00161                         continue;
00162                 //TODO: ...particularly interacting pions?
00163                 if (theSim.daughterIds().size() > 0)
00164                         continue;
00165                 answers.push_back(index);
00166                 ++index;
00167         }
00168         return answers;
00169 }
00170 
00171 void CalibratableTest::extractCandidate(const PFCandidate& cand) {
00172         if (debug_ > 3)
00173                 std::cout << "\tCandidate: "<< cand << "\n";
00174 
00175         //There may be several PFCandiates per sim particle. So we create a mini-class
00176         //to represent each one. CandidateWrapper is defined in Calibratable.
00177         //It's very easy to use, as we shall see...
00178         CandidateWrapper cw;
00179         cw.energy_ = cand.energy();
00180         cw.eta_ = cand.eta();
00181         cw.phi_ = cand.phi();
00182         cw.type_ = cand.particleId();
00183         cw.energyEcal_ = cand.ecalEnergy();
00184         cw.energyHcal_ = cand.hcalEnergy();
00185         if (debug_ > 4)
00186                 std::cout << "\t\tECAL energy = " << cand.ecalEnergy()
00187                                 << ", HCAL energy = " << cand.hcalEnergy() << "\n";
00188 
00189         //Now, extract block elements from the pfCandidate:
00190         PFCandidate::ElementsInBlocks eleInBlocks = cand.elementsInBlocks();
00191         if (debug_ > 3)
00192                 std::cout << "\tLooping over elements in blocks, "
00193                                 << eleInBlocks.size() << " of them."<< std::endl;
00194         for (PFCandidate::ElementsInBlocks::iterator bit = eleInBlocks.begin(); bit
00195                         != eleInBlocks.end(); ++bit) {
00196 
00197                 /* 
00198                  * Find PFClusters associated with this candidate.
00199                  */
00200                 
00201                 //Extract block reference
00202                 PFBlockRef blockRef((*bit).first);
00203                 //Extract index
00204                 unsigned indexInBlock((*bit).second);
00205                 //Dereference the block (what a palava!)
00206                 const PFBlock& block = *blockRef;
00207                 //And finally get a handle on the elements
00208                 const edm::OwnVector<reco::PFBlockElement> & elements = block.elements();
00209                 //get references to the candidate's track, ecal clusters and hcal clusters
00210                 switch (elements[indexInBlock].type()) {
00211                 case PFBlockElement::ECAL: {
00212                         reco::PFClusterRef clusterRef = elements[indexInBlock].clusterRef();
00213                         const PFCluster theRealCluster = *clusterRef;
00214                         CalibratableElement d(theRealCluster.energy(),
00215                                         theRealCluster.positionREP().eta(), theRealCluster.positionREP().phi(), theRealCluster.layer() );
00216                         calib_->cluster_ecal_.push_back(d);
00217                         if (debug_ > 4)
00218                                 std::cout << "\t\tECAL cluster: "<< theRealCluster << "\n";
00219 
00220                         break;
00221                 }
00222 
00223                 case PFBlockElement::HCAL: {
00224                         reco::PFClusterRef clusterRef = elements[indexInBlock].clusterRef();
00225                         const PFCluster theRealCluster = *clusterRef;
00226                         CalibratableElement d(theRealCluster.energy(),
00227                                         theRealCluster.positionREP().eta(), theRealCluster.positionREP().phi(), theRealCluster.layer() );
00228                         calib_->cluster_hcal_.push_back(d);
00229                         if (debug_ > 4)
00230                                 std::cout << "\t\tHCAL cluster: "<< theRealCluster << "\n";
00231 
00232                         break;
00233                 }
00234 
00235                 default:
00236                         if (debug_ > 3)
00237                                 std::cout << "\t\tOther block type: "
00238                                                 << elements[indexInBlock].type() << "\n";
00239                         break;
00240                 }
00241 
00242         }
00243         //For each candidate found,
00244         calib_->cands_.push_back(cw);
00245 }
00246 
00247 std::vector<unsigned> CalibratableTest::findCandidatesInDeltaR(
00248                 const PFSimParticle& pft, const std::vector<PFCandidate>& cands,
00249                 const double& deltaRCut) {
00250 
00251         unsigned index(0);
00252         std::vector<unsigned> answers;
00253 
00254         double trEta = pft.trajectoryPoint(PFTrajectoryPoint::ECALEntrance).positionREP().Eta();
00255         double trPhi = pft.trajectoryPoint(PFTrajectoryPoint::ECALEntrance).positionREP().Phi();
00256 
00257         for (std::vector<PFCandidate>::const_iterator cit = cands.begin(); cit
00258                         != cands.end(); ++cit) {
00259 
00260                 PFCandidate cand = *cit;
00261                 double cEta = cand.eta();
00262                 double cPhi = cand.phi();
00263 
00264                 if (deltaR(cEta, trEta, cPhi, trPhi) < deltaRCut) {
00265                         //accept
00266                         answers.push_back(index);
00267                 }
00268 
00269                 ++index;
00270         }
00271         return answers;
00272 }
00273 
00274 void CalibratableTest::fillTreeAndReset() {
00275         if (thisEventPasses_ && thisParticlePasses_) {
00276                 ++nParticleWrites_;
00277                 calib_->recompute();
00278                 if (debug_> 4) {
00279                         //print a summary
00280                         std::cout << *calib_;
00281                 }
00282                 tree_->Fill();
00283         } else {
00284                 ++nParticleFails_;
00285         }
00286         if (debug_ > 1)
00287                 std::cout << "\t**Finished particle.**\n";
00288         calib_->reset();
00289 }
00290 
00291 void CalibratableTest::endJob() {
00292 
00293         if (debug_> 0) {
00294                 std::cout << __PRETTY_FUNCTION__ << std::endl;
00295 
00296                 std::cout << "\tnParticleWrites: "<< nParticleWrites_
00297                                 << ", nParticleFails: "<< nParticleFails_ << "\n";
00298                 std::cout << "\tnEventWrites: "<< nEventWrites_ << ", nEventFails: "
00299                                 << nEventFails_ << "\n";
00300                 std::cout << "Leaving "<< __PRETTY_FUNCTION__ << "\n";
00301         }
00302 
00303 }
00304 
00305 double CalibratableTest::deltaR(const double& eta1, const double& eta2,
00306                 const double& phi1, const double& phi2) {
00307         double deltaEta = fabs(eta1 - eta2);
00308         double deltaPhi = fabs(phi1 - phi2);
00309         if (deltaPhi > M_PI) {
00310                 deltaPhi = 2 * M_PI- deltaPhi;
00311         }
00312         return sqrt(pow(deltaEta, 2) + pow(deltaPhi, 2));
00313 }
00314 
00315