Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018 #include <boost/ptr_container/ptr_vector.hpp>
00019 #include <boost/foreach.hpp>
00020
00021 #include <algorithm>
00022 #include <functional>
00023
00024 #include "FWCore/Framework/interface/EDProducer.h"
00025 #include "FWCore/Framework/interface/EventSetup.h"
00026 #include "FWCore/Framework/interface/ESHandle.h"
00027 #include "FWCore/Framework/interface/Event.h"
00028 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00029
00030 #include "RecoTauTag/RecoTau/interface/RecoTauBuilderPlugins.h"
00031 #include "RecoTauTag/RecoTau/interface/RecoTauCommonUtilities.h"
00032
00033 #include "DataFormats/JetReco/interface/PFJetCollection.h"
00034 #include "DataFormats/TauReco/interface/JetPiZeroAssociation.h"
00035 #include "DataFormats/TauReco/interface/PFTau.h"
00036 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00037 #include "DataFormats/Common/interface/Association.h"
00038
00039 #include "CommonTools/Utils/interface/StringCutObjectSelector.h"
00040
00041 class RecoTauProducer : public edm::EDProducer {
00042 public:
00043 typedef reco::tau::RecoTauBuilderPlugin Builder;
00044 typedef reco::tau::RecoTauModifierPlugin Modifier;
00045 typedef boost::ptr_vector<Builder> BuilderList;
00046 typedef boost::ptr_vector<Modifier> ModifierList;
00047
00048 explicit RecoTauProducer(const edm::ParameterSet& pset);
00049 ~RecoTauProducer() {}
00050 void produce(edm::Event& evt, const edm::EventSetup& es);
00051
00052 private:
00053 edm::InputTag jetSrc_;
00054 edm::InputTag jetRegionSrc_;
00055 edm::InputTag piZeroSrc_;
00056 BuilderList builders_;
00057 ModifierList modifiers_;
00058
00059 std::auto_ptr<StringCutObjectSelector<reco::PFTau> > outputSelector_;
00060
00061
00062
00063 bool buildNullTaus_;
00064 };
00065
00066 RecoTauProducer::RecoTauProducer(const edm::ParameterSet& pset) {
00067 jetSrc_ = pset.getParameter<edm::InputTag>("jetSrc");
00068 jetRegionSrc_ = pset.getParameter<edm::InputTag>("jetRegionSrc");
00069 piZeroSrc_ = pset.getParameter<edm::InputTag>("piZeroSrc");
00070
00071 typedef std::vector<edm::ParameterSet> VPSet;
00072
00073 const VPSet& builders = pset.getParameter<VPSet>("builders");
00074 for (VPSet::const_iterator builderPSet = builders.begin();
00075 builderPSet != builders.end(); ++builderPSet) {
00076
00077 const std::string& pluginType =
00078 builderPSet->getParameter<std::string>("plugin");
00079
00080 builders_.push_back(
00081 RecoTauBuilderPluginFactory::get()->create(pluginType, *builderPSet));
00082 }
00083
00084 const VPSet& modfiers = pset.getParameter<VPSet>("modifiers");
00085 for (VPSet::const_iterator modfierPSet = modfiers.begin();
00086 modfierPSet != modfiers.end(); ++modfierPSet) {
00087
00088 const std::string& pluginType =
00089 modfierPSet->getParameter<std::string>("plugin");
00090
00091 modifiers_.push_back(
00092 RecoTauModifierPluginFactory::get()->create(pluginType, *modfierPSet));
00093 }
00094
00095
00096 if (pset.exists("outputSelection")) {
00097 std::string selection = pset.getParameter<std::string>("outputSelection");
00098 if (selection != "") {
00099 outputSelector_.reset(
00100 new StringCutObjectSelector<reco::PFTau>(selection));
00101 }
00102 }
00103 buildNullTaus_ = pset.getParameter<bool>("buildNullTaus");
00104 produces<reco::PFTauCollection>();
00105 }
00106
00107 void RecoTauProducer::produce(edm::Event& evt, const edm::EventSetup& es) {
00108
00109 edm::Handle<reco::CandidateView> jetView;
00110 evt.getByLabel(jetSrc_, jetView);
00111
00112
00113 reco::PFJetRefVector jets =
00114 reco::tau::castView<reco::PFJetRefVector>(jetView);
00115
00116
00117 edm::Handle<edm::Association<reco::PFJetCollection> > jetRegionHandle;
00118 evt.getByLabel(jetRegionSrc_, jetRegionHandle);
00119
00120
00121 edm::Handle<reco::JetPiZeroAssociation> piZeroAssoc;
00122 evt.getByLabel(piZeroSrc_, piZeroAssoc);
00123
00124
00125 for (BuilderList::iterator builder = builders_.begin();
00126 builder != builders_.end(); ++builder) {
00127 builder->setup(evt, es);
00128 }
00129 for (ModifierList::iterator modifier = modifiers_.begin();
00130 modifier != modifiers_.end(); ++modifier) {
00131 modifier->setup(evt, es);
00132 }
00133
00134
00135 std::auto_ptr<reco::PFTauCollection> output(new reco::PFTauCollection());
00136
00137
00138 BOOST_FOREACH(reco::PFJetRef jetRef, jets) {
00139
00140 reco::PFJetRef jetRegionRef = (*jetRegionHandle)[jetRef];
00141 if (jetRegionRef.isNull()) {
00142 throw cms::Exception("BadJetRegionRef") << "No jet region can be"
00143 << " found for the current jet: " << jetRef.id();
00144 }
00145
00146 std::vector<reco::PFCandidatePtr> jetCands = jetRef->getPFConstituents();
00147 std::vector<reco::PFCandidatePtr> allRegionalCands =
00148 jetRegionRef->getPFConstituents();
00149
00150 std::sort(jetCands.begin(), jetCands.end());
00151 std::sort(allRegionalCands.begin(), allRegionalCands.end());
00152
00153 std::vector<reco::PFCandidatePtr> uniqueRegionalCands;
00154
00155
00156
00157 if (allRegionalCands.size() > jetCands.size())
00158 uniqueRegionalCands.reserve(allRegionalCands.size() - jetCands.size());
00159
00160
00161 std::set_difference(allRegionalCands.begin(), allRegionalCands.end(),
00162 jetCands.begin(), jetCands.end(),
00163 std::back_inserter(uniqueRegionalCands));
00164
00165
00166 const std::vector<reco::RecoTauPiZero>& piZeros = (*piZeroAssoc)[jetRef];
00167
00168 unsigned int nTausBuilt = 0;
00169 for (BuilderList::const_iterator builder = builders_.begin();
00170 builder != builders_.end(); ++builder) {
00171
00172 reco::tau::RecoTauBuilderPlugin::output_type taus(
00173 (*builder)(jetRegionRef, piZeros, uniqueRegionalCands));
00174
00175 std::for_each(taus.begin(), taus.end(),
00176 boost::bind(&reco::PFTau::setjetRef, _1, jetRef));
00177
00178
00179 output->reserve(output->size() + taus.size());
00180
00181 if (!outputSelector_.get()) {
00182 output->insert(output->end(), taus.begin(), taus.end());
00183 nTausBuilt += taus.size();
00184 } else {
00185
00186 BOOST_FOREACH(const reco::PFTau& tau, taus) {
00187 if ((*outputSelector_)(tau)) {
00188 nTausBuilt++;
00189 output->push_back(tau);
00190 }
00191 }
00192 }
00193 }
00194
00195
00196
00197 if (!nTausBuilt && buildNullTaus_) {
00198 reco::PFTau nullTau(std::numeric_limits<int>::quiet_NaN(), jetRef->p4());
00199 nullTau.setjetRef(jetRef);
00200 output->push_back(nullTau);
00201 }
00202 }
00203
00204
00205 for (reco::PFTauCollection::iterator tau = output->begin();
00206 tau != output->end(); ++tau) {
00207 for (ModifierList::const_iterator modifier = modifiers_.begin();
00208 modifier != modifiers_.end(); ++modifier) {
00209 (*modifier)(*tau);
00210 }
00211 }
00212 evt.put(output);
00213 }
00214
00215 #include "FWCore/Framework/interface/MakerMacros.h"
00216 DEFINE_FWK_MODULE(RecoTauProducer);