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(
00082 pluginType, *builderPSet));
00083
00084 }
00085
00086 const VPSet& modfiers = pset.getParameter<VPSet>("modifiers");
00087 for (VPSet::const_iterator modfierPSet = modfiers.begin();
00088 modfierPSet != modfiers.end(); ++modfierPSet) {
00089
00090 const std::string& pluginType =
00091 modfierPSet->getParameter<std::string>("plugin");
00092
00093 modifiers_.push_back(
00094 RecoTauModifierPluginFactory::get()->create(
00095 pluginType, *modfierPSet));
00096
00097 }
00098
00099
00100 if (pset.exists("outputSelection")) {
00101 std::string selection = pset.getParameter<std::string>("outputSelection");
00102 if (selection != "") {
00103 outputSelector_.reset(
00104 new StringCutObjectSelector<reco::PFTau>(selection));
00105 }
00106 }
00107 buildNullTaus_ = pset.getParameter<bool>("buildNullTaus");
00108 produces<reco::PFTauCollection>();
00109 }
00110
00111 void RecoTauProducer::produce(edm::Event& evt, const edm::EventSetup& es) {
00112
00113 edm::Handle<reco::CandidateView> jetView;
00114 evt.getByLabel(jetSrc_, jetView);
00115
00116
00117 reco::PFJetRefVector jets =
00118 reco::tau::castView<reco::PFJetRefVector>(jetView);
00119
00120
00121 edm::Handle<edm::Association<reco::PFJetCollection> > jetRegionHandle;
00122 evt.getByLabel(jetRegionSrc_, jetRegionHandle);
00123
00124
00125 edm::Handle<reco::JetPiZeroAssociation> piZeroAssoc;
00126 evt.getByLabel(piZeroSrc_, piZeroAssoc);
00127
00128
00129 for (BuilderList::iterator builder = builders_.begin();
00130 builder != builders_.end(); ++builder) {
00131 builder->setup(evt, es);
00132 }
00133 for (ModifierList::iterator modifier = modifiers_.begin();
00134 modifier != modifiers_.end(); ++modifier) {
00135 modifier->setup(evt, es);
00136 }
00137
00138
00139 std::auto_ptr<reco::PFTauCollection> output(new reco::PFTauCollection());
00140
00141
00142 BOOST_FOREACH(reco::PFJetRef jetRef, jets) {
00143
00144 reco::PFJetRef jetRegionRef = (*jetRegionHandle)[jetRef];
00145 if (jetRegionRef.isNull()) {
00146 throw cms::Exception("BadJetRegionRef") << "No jet region can be"
00147 << " found for the current jet: " << jetRef.id();
00148 }
00149
00150 std::vector<reco::PFCandidatePtr> jetCands = jetRef->getPFConstituents();
00151 std::vector<reco::PFCandidatePtr> allRegionalCands =
00152 jetRegionRef->getPFConstituents();
00153
00154 std::sort(jetCands.begin(), jetCands.end());
00155 std::sort(allRegionalCands.begin(), allRegionalCands.end());
00156
00157 std::vector<reco::PFCandidatePtr> uniqueRegionalCands;
00158
00159
00160
00161 if (allRegionalCands.size() > jetCands.size())
00162 uniqueRegionalCands.reserve(allRegionalCands.size() - jetCands.size());
00163
00164
00165 std::set_difference(allRegionalCands.begin(), allRegionalCands.end(),
00166 jetCands.begin(), jetCands.end(),
00167 std::back_inserter(uniqueRegionalCands));
00168
00169
00170 const std::vector<reco::RecoTauPiZero>& piZeros = (*piZeroAssoc)[jetRef];
00171
00172 unsigned int nTausBuilt = 0;
00173 for (BuilderList::const_iterator builder = builders_.begin();
00174 builder != builders_.end(); ++builder) {
00175
00176 reco::tau::RecoTauBuilderPlugin::output_type taus(
00177 (*builder)(jetRef, piZeros, uniqueRegionalCands));
00178
00179 std::for_each(taus.begin(), taus.end(),
00180 boost::bind(&reco::PFTau::setjetRef, _1, jetRef));
00181
00182
00183 output->reserve(output->size() + taus.size());
00184
00185 if (!outputSelector_.get()) {
00186 output->insert(output->end(), taus.begin(), taus.end());
00187 nTausBuilt += taus.size();
00188 } else {
00189
00190 BOOST_FOREACH(const reco::PFTau& tau, taus) {
00191 if ((*outputSelector_)(tau)) {
00192 nTausBuilt++;
00193 output->push_back(tau);
00194 }
00195 }
00196 }
00197 }
00198
00199
00200
00201 if (!nTausBuilt && buildNullTaus_) {
00202 reco::PFTau nullTau(std::numeric_limits<int>::quiet_NaN(), jetRef->p4());
00203 nullTau.setjetRef(jetRef);
00204 output->push_back(nullTau);
00205 }
00206 }
00207
00208
00209 for (reco::PFTauCollection::iterator tau = output->begin();
00210 tau != output->end(); ++tau) {
00211 for (ModifierList::const_iterator modifier = modifiers_.begin();
00212 modifier != modifiers_.end(); ++modifier) {
00213 (*modifier)(*tau);
00214 }
00215 }
00216 evt.put(output);
00217 }
00218
00219 #include "FWCore/Framework/interface/MakerMacros.h"
00220 DEFINE_FWK_MODULE(RecoTauProducer);