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 using namespace edm;
00021 using namespace pftools;
00022 using namespace reco;
00023
00024 CalibratableTest::CalibratableTest(const edm::ParameterSet& parameters) :
00025 debug_(0), nParticleWrites_(0), nParticleFails_(0), nEventWrites_(0), nEventFails_(0),
00026 deltaRCandToSim_(0.4) {
00027
00028 std::cout << __PRETTY_FUNCTION__ << std::endl;
00029
00030
00031 tree_ = fileservice_->make<TTree>("CalibratableTest", "");
00032 calib_ = new Calibratable();
00033 tree_->Branch("Calibratable", "pftools::Calibratable", &calib_, 32000, 2);
00034
00035 inputTagCandidates_= parameters.getParameter<InputTag>("PFCandidates");
00036 inputTagSimParticles_= parameters.getParameter<InputTag>("PFSimParticles");
00037 inputTagClustersEcal_= parameters.getParameter<InputTag>("PFClustersEcal");
00038 inputTagClustersHcal_= parameters.getParameter<InputTag>("PFClustersHcal");
00039 deltaRCandToSim_ = parameters.getParameter<double>("deltaRCandToSim");
00040 debug_= parameters.getParameter<int>("debug");
00041 }
00042
00043 CalibratableTest::~CalibratableTest() {
00044
00045 }
00046
00047 void CalibratableTest::beginJob() {
00048 if (debug_ > 1)
00049 std::cout << __PRETTY_FUNCTION__ << "\n";
00050
00051 }
00052
00053 void CalibratableTest::analyze(const edm::Event& event,
00054 const edm::EventSetup& iSetup) {
00055 if (debug_ > 1)
00056 std::cout << __PRETTY_FUNCTION__ << "\n";
00057
00058
00059 pfCandidates_ = new Handle<PFCandidateCollection>;
00060 simParticles_ = new Handle<PFSimParticleCollection>;
00061 clustersEcal_ = new Handle<PFClusterCollection>;
00062 clustersHcal_ = new Handle<PFClusterCollection>;
00063
00064 getCollection(*pfCandidates_, inputTagCandidates_, event);
00065 getCollection(*simParticles_, inputTagSimParticles_, event);
00066 getCollection(*clustersEcal_, inputTagClustersEcal_, event);
00067 getCollection(*clustersHcal_, inputTagClustersHcal_, event);
00068
00069
00070 thisEventPasses_ = true;
00071 thisParticlePasses_ = true;
00072 calib_->reset();
00073
00074 if (debug_ > 1)
00075 std::cout << "\tStarting event..."<< std::endl;
00076
00077 PFSimParticleCollection sims = **simParticles_;
00078 PFCandidateCollection candidates = **pfCandidates_;
00079 PFClusterCollection clustersEcal = **clustersEcal_;
00080 PFClusterCollection clustersHcal = **clustersHcal_;
00081
00082 if (sims.size() == 0) {
00083 std::cout << "\tNo sim particles found!" << std::endl;
00084 thisEventPasses_ = false;
00085 }
00086
00087
00088 std::vector<unsigned> primarySims = findPrimarySimParticles(sims);
00089 if (debug_) {
00090 std::cout << "\tFound "<< primarySims.size()
00091 << " primary sim particles, "<< (**pfCandidates_).size() << " pfCandidates\n";
00092 }
00093 for (std::vector<unsigned>::const_iterator cit = primarySims.begin(); cit
00094 != primarySims.end(); ++cit) {
00095
00096 if (debug_ > 1)
00097 std::cout << "\t**Starting particle...**\n";
00098 const PFSimParticle& sim = sims[*cit];
00099
00100 calib_->sim_numEvent_ = 1;
00101 calib_->sim_isMC_ = true;
00102 calib_->sim_energyEvent_ = sim.trajectoryPoint(PFTrajectoryPoint::ClosestApproach).momentum().E();
00103 calib_->sim_phi_ = sim.trajectoryPoint(PFTrajectoryPoint::ClosestApproach).momentum().Phi();
00104 calib_->sim_eta_ = sim.trajectoryPoint(PFTrajectoryPoint::ClosestApproach).momentum().Eta();
00105
00106 if (sim.nTrajectoryPoints() > PFTrajectoryPoint::ECALEntrance) {
00107 calib_->sim_etaEcal_ = sim.trajectoryPoint(PFTrajectoryPoint::ECALEntrance).positionREP().Eta();
00108 calib_->sim_phiEcal_ = sim.trajectoryPoint(PFTrajectoryPoint::ECALEntrance).positionREP().Phi();
00109 }
00110 if (sim.nTrajectoryPoints() > PFTrajectoryPoint::HCALEntrance) {
00111 calib_->sim_etaHcal_ = sim.trajectoryPoint(PFTrajectoryPoint::HCALEntrance).positionREP().Eta();
00112 calib_->sim_phiHcal_ = sim.trajectoryPoint(PFTrajectoryPoint::HCALEntrance).positionREP().Phi();
00113 }
00114
00115
00116 std::vector<unsigned> matchingCands = findCandidatesInDeltaR(sim,
00117 candidates, deltaRCandToSim_);
00118 if (debug_ > 3)
00119 std::cout << "\t\tFound candidates near sim, found "
00120 << matchingCands.size()<< " of them.\n";
00121 if (matchingCands.size() == 0)
00122 thisParticlePasses_ = false;
00123 for (std::vector<unsigned>::const_iterator mcIt = matchingCands.begin(); mcIt
00124 != matchingCands.end(); ++mcIt) {
00125 const PFCandidate& theCand = candidates[*mcIt];
00126 extractCandidate(theCand);
00127 }
00128
00129 fillTreeAndReset();
00130
00131 }
00132
00133 delete pfCandidates_;
00134 delete simParticles_;
00135 delete clustersEcal_;
00136 delete clustersHcal_;
00137
00138 if (thisEventPasses_)
00139 ++nEventWrites_;
00140 else
00141 ++nEventFails_;
00142
00143 }
00144
00145 std::vector<unsigned> CalibratableTest::findPrimarySimParticles(
00146 const std::vector<PFSimParticle>& sims) {
00147 std::vector<unsigned> answers;
00148 unsigned index(0);
00149 for (std::vector<PFSimParticle>::const_iterator cit = sims.begin(); cit
00150 != sims.end(); ++cit) {
00151 PFSimParticle theSim = *cit;
00152
00153 if (theSim.motherId() >= 0)
00154 continue;
00155 int particleId = abs(theSim.pdgCode());
00156 if (particleId != 211)
00157 continue;
00158
00159 if (theSim.daughterIds().size() > 0)
00160 continue;
00161 answers.push_back(index);
00162 ++index;
00163 }
00164 return answers;
00165 }
00166
00167 void CalibratableTest::extractCandidate(const PFCandidate& cand) {
00168 if (debug_ > 3)
00169 std::cout << "\tCandidate: "<< cand << "\n";
00170
00171
00172
00173
00174 CandidateWrapper cw;
00175 cw.energy_ = cand.energy();
00176 cw.eta_ = cand.eta();
00177 cw.phi_ = cand.phi();
00178 cw.type_ = cand.particleId();
00179 cw.energyEcal_ = cand.ecalEnergy();
00180 cw.energyHcal_ = cand.hcalEnergy();
00181 if (debug_ > 4)
00182 std::cout << "\t\tECAL energy = " << cand.ecalEnergy()
00183 << ", HCAL energy = " << cand.hcalEnergy() << "\n";
00184
00185
00186 PFCandidate::ElementsInBlocks eleInBlocks = cand.elementsInBlocks();
00187 if (debug_ > 3)
00188 std::cout << "\tLooping over elements in blocks, "
00189 << eleInBlocks.size() << " of them."<< std::endl;
00190 for (PFCandidate::ElementsInBlocks::iterator bit = eleInBlocks.begin(); bit
00191 != eleInBlocks.end(); ++bit) {
00192
00193
00194
00195
00196
00197
00198 PFBlockRef blockRef((*bit).first);
00199
00200 unsigned indexInBlock((*bit).second);
00201
00202 const PFBlock& block = *blockRef;
00203
00204 const edm::OwnVector<reco::PFBlockElement> & elements = block.elements();
00205
00206 switch (elements[indexInBlock].type()) {
00207 case PFBlockElement::ECAL: {
00208 reco::PFClusterRef clusterRef = elements[indexInBlock].clusterRef();
00209 const PFCluster theRealCluster = *clusterRef;
00210 CalibratableElement d(theRealCluster.energy(),
00211 theRealCluster.positionREP().eta(), theRealCluster.positionREP().phi(), theRealCluster.layer() );
00212 calib_->cluster_ecal_.push_back(d);
00213 if (debug_ > 4)
00214 std::cout << "\t\tECAL cluster: "<< theRealCluster << "\n";
00215
00216 break;
00217 }
00218
00219 case PFBlockElement::HCAL: {
00220 reco::PFClusterRef clusterRef = elements[indexInBlock].clusterRef();
00221 const PFCluster theRealCluster = *clusterRef;
00222 CalibratableElement d(theRealCluster.energy(),
00223 theRealCluster.positionREP().eta(), theRealCluster.positionREP().phi(), theRealCluster.layer() );
00224 calib_->cluster_hcal_.push_back(d);
00225 if (debug_ > 4)
00226 std::cout << "\t\tHCAL cluster: "<< theRealCluster << "\n";
00227
00228 break;
00229 }
00230
00231 default:
00232 if (debug_ > 3)
00233 std::cout << "\t\tOther block type: "
00234 << elements[indexInBlock].type() << "\n";
00235 break;
00236 }
00237
00238 }
00239
00240 calib_->cands_.push_back(cw);
00241 }
00242
00243 std::vector<unsigned> CalibratableTest::findCandidatesInDeltaR(
00244 const PFSimParticle& pft, const std::vector<PFCandidate>& cands,
00245 const double& deltaRCut) {
00246
00247 unsigned index(0);
00248 std::vector<unsigned> answers;
00249
00250 double trEta = pft.trajectoryPoint(PFTrajectoryPoint::ECALEntrance).positionREP().Eta();
00251 double trPhi = pft.trajectoryPoint(PFTrajectoryPoint::ECALEntrance).positionREP().Phi();
00252
00253 for (std::vector<PFCandidate>::const_iterator cit = cands.begin(); cit
00254 != cands.end(); ++cit) {
00255
00256 PFCandidate cand = *cit;
00257 double cEta = cand.eta();
00258 double cPhi = cand.phi();
00259
00260 if (deltaR(cEta, trEta, cPhi, trPhi) < deltaRCut) {
00261
00262 answers.push_back(index);
00263 }
00264
00265 ++index;
00266 }
00267 return answers;
00268 }
00269
00270 void CalibratableTest::fillTreeAndReset() {
00271 if (thisEventPasses_ && thisParticlePasses_) {
00272 ++nParticleWrites_;
00273 calib_->recompute();
00274 if (debug_> 4) {
00275
00276 std::cout << *calib_;
00277 }
00278 tree_->Fill();
00279 } else {
00280 ++nParticleFails_;
00281 }
00282 if (debug_ > 1)
00283 std::cout << "\t**Finished particle.**\n";
00284 calib_->reset();
00285 }
00286
00287 void CalibratableTest::endJob() {
00288
00289 if (debug_> 0) {
00290 std::cout << __PRETTY_FUNCTION__ << std::endl;
00291
00292 std::cout << "\tnParticleWrites: "<< nParticleWrites_
00293 << ", nParticleFails: "<< nParticleFails_ << "\n";
00294 std::cout << "\tnEventWrites: "<< nEventWrites_ << ", nEventFails: "
00295 << nEventFails_ << "\n";
00296 std::cout << "Leaving "<< __PRETTY_FUNCTION__ << "\n";
00297 }
00298
00299 }
00300
00301 double CalibratableTest::deltaR(const double& eta1, const double& eta2,
00302 const double& phi1, const double& phi2) {
00303 double deltaEta = fabs(eta1 - eta2);
00304 double deltaPhi = fabs(phi1 - phi2);
00305 if (deltaPhi > M_PI) {
00306 deltaPhi = 2 * M_PI- deltaPhi;
00307 }
00308 return sqrt(pow(deltaEta, 2) + pow(deltaPhi, 2));
00309 }
00310
00311